Arcane  v4.1.1.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
ZoltanMeshPartitioner.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2025 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/* Partitioneur de maillage utilisant la bibliotheque Zoltan. */
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// Au cas ou on utilise mpich2 ou 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: L'implémentation actuelle ne supporte que les uniqueId() sur 32 bits.
59
60/*---------------------------------------------------------------------------*/
61/*---------------------------------------------------------------------------*/
62
63namespace Arcane
64{
65
66/*---------------------------------------------------------------------------*/
67/*---------------------------------------------------------------------------*/
68
69
70enum ZoltanModel
71{
72 ZOLTAN_HG = (1<<0),
73 ZOLTAN_GRAPH = (1<<1),
74 ZOLTAN_GEOM = (1<<2)
75};
76
77/*---------------------------------------------------------------------------*/
82class ZoltanInfo
83 : public TraceAccessor
84{
85public:
86 ZoltanInfo(MeshPartitionerBase* basePartitionner,ostream* ofile, int model=1, const Real edgeWeightMultiplier = 1.)
87 : TraceAccessor(basePartitionner->mesh()->traceMng())
88 , m_mesh_partitionner_base(basePartitionner)
89 , m_nbEdges(0)
90 , m_nbPins(0)
91 , m_ofile(ofile)
92 , m_model(model)
93 , m_edgeGIDStart(0)
94 , m_edge_weight_multiplier(edgeWeightMultiplier)
95 {
96 m_own_cells = m_mesh_partitionner_base->mesh()->ownCells();
97 build();
98 }
99 MeshPartitionerBase* m_mesh_partitionner_base;
100private:
101 CellGroup m_own_cells;
102 int m_nbEdges;
103 int m_nbPins;
104 ostream* m_ofile;
105 int m_model;
106 int m_edgeGIDStart;
107 Real m_edge_weight_multiplier;
108 std::set<std::pair<Int64, Int64> > m_weight_set;
109public:
110 void build()
111 {
112 info() << "ZoltanInfo::build()";
113 if (m_model & ZOLTAN_GEOM) // No topological informations to build
114 return;
115
116 Integer nbOwnCells = m_own_cells.size();
117 Integer nbOwnEdges = 0;
118
119 if (m_model & ZOLTAN_HG) {
120 nbOwnEdges += nbOwnCells;
121 }
122 if (m_model & ZOLTAN_GRAPH) {
123 // Neighboors, each and then neighboorhood HE
124 nbOwnEdges += nbOwnCells * 6;
125 }
126
127 // Compute m_edgeGIDStart = Sum_k<i nbOwnEdges
128 IParallelMng* pm = m_mesh_partitionner_base->mesh()->parallelMng();
129 Int32UniqueArray scanouille(pm->commSize());
130 scanouille.fill(0);
131 scanouille[pm->commRank()] = nbOwnEdges;
132 pm->scan(Parallel::ReduceSum,scanouille);
133 m_edgeGIDStart = 0;
134 for (int i = 0 ; i < pm->commRank() ; ++i) {
135 m_edgeGIDStart += scanouille[i];
136 }
137
138 m_nbEdges = 0;
139 m_nbPins = 0;
140 ENUMERATE_CELL(iCell,m_own_cells)
141 {
142 if (!m_mesh_partitionner_base->cellUsedWithConstraints(*iCell))
143 continue;
144
145 int nbNgb = m_mesh_partitionner_base->nbNeighbourCellsWithConstraints(*iCell);
146 if (m_model & ZOLTAN_HG) {
147 m_nbEdges ++;
148 m_nbPins += nbNgb + 1;
149 }
150 if (m_model & ZOLTAN_GRAPH) {
151 m_nbEdges += nbNgb;
152 m_nbPins += (nbNgb*2);
153 }
154 }
155
156 info() << "nbEdges=" << m_nbEdges << " ; nbPins=" << m_nbPins;
157
158 if(m_mesh_partitionner_base->haveWeakConstraints())
159 {
160 MeshVariableScalarRefT<Face, Integer> weak_constraint(VariableBuildInfo(m_mesh_partitionner_base->mesh(), "EdgeWeight"));
161 ENUMERATE_FACE(iface,m_mesh_partitionner_base->mesh()->ownFaces())
162 {
163 const Face& face = *iface;
164 const Cell& bCell = iface->backCell();
165 const Cell& fCell = iface->frontCell();
166 if(bCell.null() || fCell.null())
167 continue;
168 if(weak_constraint[face]==2)
169 {
170 m_weight_set.insert(std::pair<Int64,Int64>(face.backCell().uniqueId(), face.frontCell().uniqueId()));
171 m_weight_set.insert(std::pair<Int64,Int64>(face.frontCell().uniqueId(), face.backCell().uniqueId()));
172 }
173 }
174 }
175 } // end build()
176
177
178public:
179 static int getHgNumVertices(void *data, int *ierr)
180 {
181 ARCANE_UNUSED(ierr);
182
183 ZoltanInfo* zi = (ZoltanInfo*)data;
184
185 /*
186 * Supply this query function to Zoltan with Zoltan_Set_Num_Obj_Fn().
187 * It returns the number of vertices that this process owns.
188 *
189 * The parallel hypergraph method requires that vertex global IDs and
190 * weights are returned by the application with query functions. The
191 * global IDs should be unique, and no two processes should
192 * return the same vertex.
193 */
194 return zi->m_mesh_partitionner_base->nbOwnCellsWithConstraints();
195 }
196
197 static void getHgVerticesAndWeights(void *data, int num_gid_entries,
198 int num_lid_entries, ZOLTAN_ID_PTR gids, ZOLTAN_ID_PTR lids,
199 int wgt_dim, float *obj_weights, int *ierr)
200 {
201 ARCANE_UNUSED(num_gid_entries);
202 ARCANE_UNUSED(num_lid_entries);
203
204 ZoltanInfo* zi = (ZoltanInfo*)data;
205 /*
206 * Supply this query function to Zoltan with Zoltan_Set_Obj_List_Fn().
207 *
208 * It supplies vertex global IDs, local IDs,
209 * and weights to the Zoltan library. The application has previously
210 * indicated the number of weights per vertex by setting the
211 * OBJ_WEIGHT_DIM parameter.
212 */
213 //int nb_cell = zi->m_own_cells.size();
214
215 SharedArray<float> cells_weights;
216 if (wgt_dim != 0) { // No weight, does not mean no limitation like for Arcane
217 if (zi->m_mesh_partitionner_base->nbCellWeight() <= wgt_dim) { // We can use all criteria, so we do
218 cells_weights = zi->m_mesh_partitionner_base->cellsWeightsWithConstraints(wgt_dim);
219 }
220 else { // We need more criteria than available
221 // So we try to balance memory !
222 cells_weights = zi->m_mesh_partitionner_base->cellsSizeWithConstraints();
223 }
224 ArrayView<float> view_weights(cells_weights.size(), obj_weights);
225 view_weights.copy(cells_weights);
226 }
227
228 int index = 0;
229 ENUMERATE_CELL(icell,zi->m_own_cells)
230 {
231 if (!zi->m_mesh_partitionner_base->cellUsedWithConstraints(*icell))
232 continue;
233
234 gids[index] = (*icell).uniqueId().asInt32();
235 lids[index] = zi->m_mesh_partitionner_base->localIdWithConstraints(*icell);
236
237 if (zi->m_ofile) {
238 for( int j=0; j< wgt_dim; ++j ){
239 float weight = cells_weights[lids[index]*wgt_dim+j];
240 *obj_weights++ = weight;
241 // zi->info() << " Weight uid=" << gids[index] << " w=" << weight;
242 (*zi->m_ofile) << " Weight uid=" << gids[index] << " w=" << weight << '\n';
243 }
244 }
245 index ++;
246 }
247 *ierr = ZOLTAN_OK;
248 }
249
250 static void getHgSizeAndFormat(void *data, int *num_lists, int *num_pins, int *format, int *ierr)
251 {
252 ZoltanInfo* zi = (ZoltanInfo*)data;
253
254 zi->info() << " ZoltanInfo::getHgSizeAndFormat() ";
255
256 /*
257 * Supply this query function to Zoltan with Zoltan_Set_HG_Size_CS_Fn().
258 * It tells Zoltan the number of rows or columns to be supplied by
259 * the process, the number of pins (non-zeroes) in the rows or columns,
260 * and whether the pins are provided in compressed row format or
261 * compressed column format.
262 */
263 *format = ZOLTAN_COMPRESSED_EDGE;
264 *num_pins = zi->m_nbPins;
265 *num_lists = zi->m_nbEdges;
266
267 // if (zi->m_model == ZOLTAN_MY_HG) {
268 // *format = ZOLTAN_COMPRESSED_VERTEX;
269 // }
270 zi->info() << " ZoltanInfo::getHgSizeAndFormat() " << " num_list= " << *num_lists << " num_pins= " << *num_pins;
271 *ierr =ZOLTAN_OK;
272 }
273
274 static void getHg(void *data, int num_gid_entries,
275 int nrowcol, int npins, int format,
276 ZOLTAN_ID_PTR z_vtxedge_GID, int *z_vtxedge_ptr, ZOLTAN_ID_PTR z_pin_GID, int *ierr)
277 {
278 ZoltanInfo* zi = (ZoltanInfo*)data;
279
280 zi->info() << " ZoltanInfo::getHg() "
281 << " num_gid_entries= " << num_gid_entries
282 << " nrowcol= " << nrowcol
283 << " npins= " << npins
284 << " format= " << format;
285
286 // 6 neighbors max ? >> Not true for groups !
287 Int64UniqueArray neighbour_cells;
288 neighbour_cells.reserve(6);
289 int indexPin = 0;
290 int indexEdge = 0;
291 int gid = zi->m_edgeGIDStart;
292 z_vtxedge_ptr[0] = 0;
293 ENUMERATE_CELL(i_item, zi->m_own_cells) {
294 const Cell& item = *i_item;
295
296 if (!zi->m_mesh_partitionner_base->cellUsedWithConstraints(item))
297 continue;
298
299 neighbour_cells.resize(0);
300 zi->m_mesh_partitionner_base->getNeighbourCellsUidWithConstraints(item,
301 neighbour_cells);
302
303 if (zi->m_model & ZOLTAN_HG)
304 {
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 {
316 // size 2 edges are face communications
317 // other are neighborhood description
318 for( Integer z=0; z<neighbour_cells.size(); ++z ) {
319 z_pin_GID[indexPin++] = CheckedConvert::toInteger(neighbour_cells[z]);
320 z_pin_GID[indexPin++] = (item.uniqueId().asInt32());
321 z_vtxedge_GID[indexEdge] = gid++;
322 z_vtxedge_ptr[indexEdge+1] = z_vtxedge_ptr[indexEdge] + 2;
323 indexEdge++;
324 }
325 }
326 }
327 if (zi->m_ofile) {
328 for (int i = 0 ; i < nrowcol ; ++i) {
329 (*zi->m_ofile) << "*topo GID " << z_vtxedge_GID[i]
330 << " : " << z_vtxedge_ptr[i+1] - z_vtxedge_ptr[i];
331 for (int j = z_vtxedge_ptr[i] ; j < z_vtxedge_ptr[i+1] ; ++j) {
332 (*zi->m_ofile) << '\t' << z_pin_GID[j];
333 }
334 (*zi->m_ofile) << '\n';
335 }
336 }
337 *ierr = ZOLTAN_OK;
338 }
339
340
341 /*
342 For a list of objects, it returns the per-objects sizes (in bytes)
343 of the data buffers needed to pack object data.
344 void (*)(void*, int, int, int, ZOLTAN_ID_TYPE*, ZOLTAN_ID_TYPE*, int*, int*)
345 */
346
347 static void getHgVertexSizes (void *data, int num_gid_entries, int num_lid_entries,
348 int num_ids, ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids,
349 int *sizes, int *ierr)
350 {
351 ARCANE_UNUSED(num_gid_entries);
352 ARCANE_UNUSED(num_lid_entries);
353 ARCANE_UNUSED(global_ids);
354
355 ZoltanInfo* zi = (ZoltanInfo*)data;
356
357 SharedArray<float> cellSizes = zi->m_mesh_partitionner_base->cellsSizeWithConstraints();
358
359 for (int i = 0 ; i < num_ids ; ++i)
360 {
361 sizes[i] = Convert::toInteger(cellSizes[local_ids[i]]);
362 }
363
364 *ierr = ZOLTAN_OK;
365 }
366
367 static void getHgEdgeWeightSize(void *data, int *num_edges, int *ierr)
368 {
369 ZoltanInfo* zi = (ZoltanInfo*)data;
370
371 /*
372 * Supply this query function to Zoltan with Zoltan_Set_HG_Size_Edge_Weights_Fn().
373 * It tells Zoltan the number edges for which this process will supply
374 * edge weights. The number of weights per edge was specified with the
375 * parameter EDGE_WEIGHT_DIM. If more than one process provides a weight
376 * for the same edge, the multiple weights are resolved according the
377 * value for the PHG_EDGE_WEIGHT_OPERATION parameter.
378 */
379
380 *num_edges = zi->m_nbEdges;
381 *ierr = ZOLTAN_OK;
382 }
383
384 static void getHgEdgeWeights(void *data, int num_gid_entries,
385 int num_lid_entries, int nedges, int edge_weight_dim,
386 ZOLTAN_ID_PTR edge_GID, ZOLTAN_ID_PTR edge_LID, float *edge_weight, int *ierr)
387 {
388 ARCANE_UNUSED(num_lid_entries);
389
390 ZoltanInfo* zi = (ZoltanInfo*)data;
391 /*
392 * Supply this query function to Zoltan with Zoltan_Set_HG_Edge_Weights_Fn().
393 * It tells Zoltan the weights for some subset of edges.
394 */
395
396 zi->info() << " ZoltanInfo::getHgEdgeWeights() "
397 << " num_gid_entries= " << num_gid_entries
398 << " nedges= " << nedges
399 << " edge_weight_dim= " << edge_weight_dim;
400
401 UniqueArray<float> connectivityWeights;
402 Int64UniqueArray neighbour_cells;
403 neighbour_cells.reserve(6);
404 connectivityWeights.reserve(6);
405 int indexEdge = 0;
406
407 ENUMERATE_CELL(i_item, zi->m_own_cells) {
408 Cell item = *i_item;
409
410 if (!zi->m_mesh_partitionner_base->cellUsedWithConstraints(item))
411 continue;
412
413 connectivityWeights.resize(0);
414 neighbour_cells.resize(0);
415 bool hg_model=(zi->m_model & ZOLTAN_HG);
416 Real he_weight= 0;
417 he_weight =
418 zi->m_mesh_partitionner_base->getNeighbourCellsUidWithConstraints(item,
419 neighbour_cells, &connectivityWeights, hg_model);
420
421 if (zi->m_model & ZOLTAN_HG) {
422 if (!(zi->m_model & ZOLTAN_GRAPH)) {
423 // We have to sum-up de Weight of pins to define HyperEdge's one.
424 for( Integer z=0; z<neighbour_cells.size(); ++z ) {
425 he_weight += connectivityWeights[z];
426 }
427 }
428 edge_GID[indexEdge] = indexEdge + zi->m_edgeGIDStart;
429 edge_LID[indexEdge] = indexEdge;
430 edge_weight[indexEdge++] = (float)he_weight;
431 }
432
433 if (zi->m_model & ZOLTAN_GRAPH) {
434 // size 2 edges are face communications
435 // other are neighborhood description
436 for( Integer z=0; z<neighbour_cells.size(); ++z ){
437 edge_GID[indexEdge] = indexEdge + zi->m_edgeGIDStart;
438 edge_LID[indexEdge] = indexEdge;
439 Real w = connectivityWeights[z];
440 edge_weight[indexEdge] = (float)w;
441 if (zi->m_mesh_partitionner_base->haveWeakConstraints()){
442 std::pair<Int64,Int64> items(item.uniqueId(), neighbour_cells[z]);
443 if(zi->m_mesh_partitionner_base->cellUsedWithWeakConstraints(items)){
444 edge_weight[indexEdge] *= (float)zi->m_edge_weight_multiplier;
445 }
446 }
447 ++indexEdge;
448 }
449 }
450 }
451 *ierr = ZOLTAN_OK;
452 }
453
454 // Return the dimension of a vertex, for geometric methods
455 static int get_num_geometry(void *data, int *ierr)
456 {
457 ARCANE_UNUSED(data);
458 *ierr = ZOLTAN_OK;
459 return 3;
460 }
461
462 // Return the coordinates of my objects (vertices), for geometric methods.
463 static void get_geometry_list(void *data, int sizeGID, int sizeLID,
464 int num_obj,
465 ZOLTAN_ID_PTR global_ids, ZOLTAN_ID_PTR local_ids,
466 int num_dim, double *geom_vec, int *ierr)
467 {
468 ARCANE_UNUSED(global_ids);
469 ARCANE_UNUSED(local_ids);
470
471 ZoltanInfo* zi = (ZoltanInfo*)data;
472
473 if ( (sizeGID != 1) || (sizeLID != 1) || (num_dim > 3)){
474 *ierr = ZOLTAN_FATAL;
475 return;
476 }
477
478 VariableNodeReal3& coords(zi->m_mesh_partitionner_base->mesh()->nodesCoordinates());
479
480
481 int i=0;
482 ENUMERATE_CELL(icell,zi->m_own_cells){
483 if (!zi->m_mesh_partitionner_base->cellUsedWithConstraints(*icell))
484 continue;
485
486 // on calcul un barycentre
487
488 Real3 bar;
489
490 for( Integer z=0, zs = (*icell).nbNode(); z<zs; ++z ){
491 const Node& node = (*icell).node(z);
492 bar += coords[node];
493 }
494 bar /= Convert::toDouble((*icell).nbNode());
495
496 geom_vec[num_dim*i ] = bar.x;
497 geom_vec[num_dim*i+1] = bar.y;
498 geom_vec[num_dim*i+2] = bar.z;
499 i += 1;
500 }
501
502 String s = platform::getEnvironmentVariable("ZOLTAN_MODEL");
503 if (!s.null() && (s == "MYGEOM")) {
504 for (int i = 0 ; i < num_obj ; ++i) {
505 geom_vec[num_dim*i+num_dim-1] = 0.0;
506 }
507 }
508
509 *ierr = ZOLTAN_OK;
510 return;
511 }
512};
513
514/*---------------------------------------------------------------------------*/
515/*---------------------------------------------------------------------------*/
519class ZoltanMeshPartitioner
521{
522 public:
523
524 explicit ZoltanMeshPartitioner(const ServiceBuildInfo& sbi);
525
526 public:
527
528 virtual void build() {}
529
530 public:
531
532 virtual void partitionMesh(bool initial_partition);
533 virtual void partitionMesh(bool initial_partition,Int32 nb_part);
534
535 virtual void notifyEndPartition();
536
537 private:
538};
539
540/*---------------------------------------------------------------------------*/
541/*---------------------------------------------------------------------------*/
542
543ZoltanMeshPartitioner::
544ZoltanMeshPartitioner(const ServiceBuildInfo& sbi)
545 : ArcaneZoltanMeshPartitionerObject(sbi)
546{
547}
548
549/*---------------------------------------------------------------------------*/
550/*---------------------------------------------------------------------------*/
551
553partitionMesh(bool initial_partition)
554{
555 Int32 nb_part = mesh()->parallelMng()->commSize();
556 partitionMesh(initial_partition,nb_part);
557}
558
559/*---------------------------------------------------------------------------*/
560/*---------------------------------------------------------------------------*/
561
563partitionMesh(bool initial_partition,Int32 nb_part)
564{
565 info() << "Load balancing with Zoltan\n";
566
567 float ver;
568
569 int rc = ::Zoltan_Initialize(0,0,&ver);
570 Zoltan_Memory_Debug(2);
571 if (rc != ZOLTAN_OK)
572 fatal() << "Can not initialize zoltan (r=" << rc << ")";
573 IMesh* mesh = this->mesh();
574 IParallelMng* pm = mesh->parallelMng();
575 Integer nb_rank = pm->commSize();
576
577 if (nb_part<nb_rank)
578 throw ArgumentException(A_FUNCINFO,"partition with nb_part<nb_rank");
579
580 if (nb_rank==1){
581 warning() << "Unable to test load balancing on a single sub-domain";
582 return;
583 }
584
585 Integer nb_weight = nbCellWeight();
586
587 struct Zoltan_Struct *zz;
588 int changes;
589 int numGidEntries;
590 int numLidEntries;
591 int numImport;
592 ZOLTAN_ID_PTR importGlobalIds;
593 ZOLTAN_ID_PTR importLocalIds;
594 int *importProcs;
595 int *importToPart;
596 int numExport;
597 ZOLTAN_ID_PTR exportGlobalIds;
598 int *exportProcs;
599
600 /******************************************************************
601 ** Prepare to partition the example hypergraph using the Zoltan
602 ** parallel hypergraph package.
603 **
604 ** Set parameter values, and supply to Zoltan the query functions
605 ** defined in the example library which will provide to the Zoltan
606 ** library the hypergraph data.
607 ******************************************************************/
608 zz = Zoltan_Create(*(MPI_Comm*)getCommunicator());
609
610 Zoltan_Set_Param(zz, "RCB_REUSE", "1"); // Try to do "repartitioning"
611 // option entre PHG et RCB (hypergraphe et coordonnee bissection)
612
613 // Anciennes valeurs par defaut si pas de variable d'environnement (cf OLD_HG)
614 bool usePHG = true;
615
616 String s = platform::getEnvironmentVariable("ZOLTAN_MODEL");
617
618 if (options()) {
619 if (!s.null() && options()->model.isPresent())
620 fatal() << "Conflicting configuration between ZOLTAN_MODEL environment variable and user data set";
621 usePHG = options()->useHypergraphe();
622 if (!usePHG)
623 s = "RCB";
624 else if (s.null())
625 s = options()->model();
626 } else if (s.null()) {
627 s = "OLDHG";
628 }
629 // => s != null
630
631 if (s == "HYBRID") {
632 if (subDomain()->commonVariables().globalIteration() <= 2)
633 s = "RCB";
634 else
635 s = "DUALHG";
636 }
637
638 String algo = "HYPERGRAPH";
639 if (s == "MYHG" || s == "OLDHG" /* nuance avec OLD_HG ? */ || s == "DUALHG" || s == "GRAPH") {
640 algo = "HYPERGRAPH";
641 usePHG = true;
642 } else if (s == "RIB" || s== "HSFC") {
643 algo = s;
644 usePHG = false;
645 } else if (s == "RCB" || s == "MYGEOM") {
646 algo = "RCB";
647 usePHG = false;
648 } else {
649 fatal() << "Undefined zoltan model '" << s << "'";
650 }
651
652 int model = ZOLTAN_HG;
653 if (usePHG == false) {
654 model = ZOLTAN_GEOM;
655 }
656
657 if (s == "DUALHG") {
658 model |= ZOLTAN_GRAPH;
659 } else if (s == "GRAPH") {
660 model = ZOLTAN_GRAPH;
661 }
662
663 if (usePHG) {
664 nb_weight = 1; // up to 2 weights for RCB
665 }
666
667 Real edgeWeightMultiplier = 1;
668 Integer repartitionFrequency = 10;
669 Real imbalanceTol = 1.05;
670 Real phgRepartMultiplier = 10;
671 Integer phgOutputLevel = 0;
672 Integer debugLevel = 0;
673
674 if(options())
675 {
676 edgeWeightMultiplier = options()->edgeWeightMultiplier();
677 repartitionFrequency = options()->repartFrequency();
678 imbalanceTol = options()->imbalanceTol();
679 phgRepartMultiplier = options()->phgRepartMultiplier();
680 phgOutputLevel = options()->phgOutputLevel();
681 debugLevel = options()->debugLevel();
682 }
683
684 /* General parameters */
685 info() << "Zoltan: utilise un repartitionnement " << algo <<" (" << s << ").";
686 Zoltan_Set_Param(zz, "LB_METHOD", algo.localstr()); /* partitioning method */
687 Zoltan_Set_Param(zz, "HYPERGRAPH_PACKAGE", "PHG"); /* version of method */
688 //if(!initial_partition)
689 if(mesh->subDomain()->commonVariables().globalIteration()==1)
690 {
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 {
696 Zoltan_Set_Param(zz, "LB_APPROACH", "REPARTITION"); /* version of method */
697 info() << "Zoltan: Repartition";
698 }
699 else
700 {
701 Zoltan_Set_Param(zz, "LB_APPROACH", "REFINE"); /* version of method */
702 info() << "Zoltan: Refine";
703 }
704
705 String s_imbalaceTol(String::fromNumber(imbalanceTol));
706 Zoltan_Set_Param(zz, "IMBALANCE_TOL", s_imbalaceTol.localstr()); /* version of method */
707
708 String s_nb_part(String::fromNumber(nb_part));
709 Zoltan_Set_Param(zz, "NUM_GLOBAL_PARTS", s_nb_part.localstr());
710
711 Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1");/* global IDs are integers */
712 Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "1");/* local IDs are integers */
713 Zoltan_Set_Param(zz, "RETURN_LISTS", "EXPORT"); /* only export lists */
714 String s_nb_weight(String::fromNumber(nb_weight));
715
716 Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", s_nb_weight.localstr()); /* 2 vtx weight */
717
718 /* Graph parameters */
719 Zoltan_Set_Param(zz, "ADD_OBJ_WEIGHT", "NONE"); /* Don't calculate extra weights */
720
721 if (usePHG) {
722 /* Parallel hypergraph parameters */
723 Zoltan_Set_Param(zz, "FINAL_OUTPUT", "0"); /* provide stats at the end */
724 Zoltan_Set_Param(zz, "PHG_USE_TIMERS", "1");
725
726 String s_phgOutputLevel(String::fromNumber(phgOutputLevel));
727 Zoltan_Set_Param(zz, "PHG_OUTPUT_LEVEL", s_phgOutputLevel.localstr());
728
729 // Utilisation en debug
730 Zoltan_Set_Param(zz, "CHECK_HYPERGRAPH", "0"); /* see User's Guide */
731 String s_debugLevel(String::fromNumber(debugLevel));
732 Zoltan_Set_Param(zz, "DEBUG_LEVEL", s_debugLevel.localstr());
733
734 // La methode par defaut (ipm) donne en theorie de meilleurs resultats
735 // mais elle semble des fois bloquer avec la version 2.0 de zoltan
736
737 Zoltan_Set_Param(zz, "PHG_COARSENING_METHOD", "IPM");
738 //Zoltan_Set_Param(zz, "PHG_COARSEPARTITION_METHOD", "GREEDY");
739 //Zoltan_Set_Param(zz, "PHG_CUT_OBJECTIVE", "HYPEREDGES");
740
741 Zoltan_Set_Param(zz, "EDGE_WEIGHT_DIM", "1");
742 Zoltan_Set_Param(zz, "PHG_EDGE_WEIGHT_OPERATION", "add");
743
744 if(!initial_partition)
745 {
746 // Number of iterations between 2 load balances. Used to scale the cost
747 // of data migration.
748 String s_phgRepartMultiplier(String::fromNumber(phgRepartMultiplier));
749 Zoltan_Set_Param(zz, "PHG_REPART_MULTIPLIER", s_phgRepartMultiplier.localstr());
750 }
751
752 }
753 else {
754 Zoltan_Set_Param(zz, "KEEP_CUTS", "0");
755 Zoltan_Set_Param(zz, "RCB_OUTPUT_LEVEL", "0");
756 Zoltan_Set_Param(zz, "RCB_RECTILINEAR_BLOCKS", "0");
757 // The 2 following options are for multi-weights partitioning
758 // Object weights are not comparable
759 Zoltan_Set_Param(zz, "OBJ_WEIGHTS_COMPARABLE", "0");
760 // 1 : minimize total (ie Sum over all phases)
761 // 2 : between
762 // 3 : Minimize worst case for each phase
763 Zoltan_Set_Param(zz, "RCB_MULTICRITERIA_NORM", "3");
764 Zoltan_Set_Param(zz, "RCB_MAX_ASPECT_RATIO", "10");
765 }
766
767 bool dump_infos = false;
768 if (platform::getEnvironmentVariable("ARCANE_DEBUG_PARTITION")=="TRUE")
769 dump_infos = true;
770 ofstream ofile;
771 if (dump_infos){
772 StringBuilder fname;
773 Integer iteration = mesh->subDomain()->commonVariables().globalIteration();
774 fname = "weigth-";
775 fname += pm->commRank();
776 fname += "-";
777 fname += iteration;
778 String f(fname);
779 ofile.open(f.localstr());
780 }
781
782 /* Application defined query functions (defined in exphg.c) */
783
784 // initialisation pour la gestion des contraintes
785 initConstraints();
786
787 ostream* zofile = 0;
788 if (dump_infos)
789 zofile = &ofile;
790
791 ScopedPtrT<ZoltanInfo> zoltan_info(new ZoltanInfo(this,zofile, model, edgeWeightMultiplier));
792 Zoltan_Set_Num_Obj_Fn(zz,&ZoltanInfo::getHgNumVertices, zoltan_info.get());
793 Zoltan_Set_Obj_List_Fn(zz, &ZoltanInfo::getHgVerticesAndWeights, zoltan_info.get());
794 if (!initial_partition) {
795 Zoltan_Set_Obj_Size_Multi_Fn(zz, &ZoltanInfo::getHgVertexSizes,
796 zoltan_info.get());
797 }
798 if (usePHG){
799 info() <<"Setting up HG Callbacks";
800 Zoltan_Set_HG_Size_CS_Fn(zz, &ZoltanInfo::getHgSizeAndFormat, zoltan_info.get());
801 Zoltan_Set_HG_CS_Fn(zz, &ZoltanInfo::getHg, zoltan_info.get());
802 Zoltan_Set_HG_Size_Edge_Weights_Fn(zz, &ZoltanInfo::getHgEdgeWeightSize, zoltan_info.get());
803 Zoltan_Set_HG_Edge_Weights_Fn(zz, &ZoltanInfo::getHgEdgeWeights, zoltan_info.get());
804 }
805 else
806 {
807 info() <<"Setting up Geom Callbacks";
808 Zoltan_Set_Num_Geom_Fn(zz, ZoltanInfo::get_num_geometry, zoltan_info.get());
809 Zoltan_Set_Geom_Multi_Fn(zz, ZoltanInfo::get_geometry_list, zoltan_info.get());
810 }
811
812 /* Parallel partitioning occurs now */
813
814 int* export_partitions = 0;
815 ZOLTAN_ID_PTR export_local_ids = 0;
816
817 info() << "Doing partition";
818 rc = Zoltan_LB_Partition(zz, &changes, &numGidEntries, &numLidEntries,
819 &numImport, &importGlobalIds, &importLocalIds,
820 &importProcs, &importToPart,
821 &numExport, &exportGlobalIds, &export_local_ids,
822 &exportProcs, &export_partitions);
823 pm->barrier();
824
825 VariableItemInt32& cells_new_owner = mesh->toPrimaryMesh()->itemsNewOwner(IK_Cell);
826 ENUMERATE_ITEM(icell,mesh->ownCells()){
827 cells_new_owner[icell] = (*icell).owner();
828 }
829
830 invertArrayLid2LidCompacted();
831
832 {
833 CellInfoListView items_internal(m_cell_family);
834 Integer nb_export = numExport;
835 info() <<"numExport = "<<numExport;
836 for( Integer i=0; i<nb_export; ++i ){
837 Item item = items_internal[ localIdWithConstraints(export_local_ids[i]) ];
838 // changement pour la maille ou tout le groupe s'il y a lieu
839 changeCellOwner(item, cells_new_owner, export_partitions[i]);
840 if (dump_infos) // ces infos ne tiennent pas compte des groupes/contraintes
841 ofile << "EXPORT: uid=" << ItemPrinter(item) << " OLD=" << item.owner()
842 << " NEW=" << cells_new_owner[item] << " PROC=" << exportProcs[i] << '\n';
843 }
844 // pinfo() << "Proc " << my_rank << " nombre de mailles changeant de domaine: "
845 // << nb_export << " changes=" << changes
846 // << " ZoltanMem=" << Zoltan_Memory_Usage(ZOLTAN_MEM_STAT_MAXIMUM);
847 }
848
849 // Libere la memoire de zoltan_info
850 zoltan_info = 0;
851 if (dump_infos)
852 ofile.close();
853
854 if (numExport < 0) {
855 sleep(3);
856 MPI_Abort(MPI_COMM_WORLD, EXIT_FAILURE);
857 }
858
859 Zoltan_LB_Free_Part(&exportGlobalIds, &export_local_ids,
860 &exportProcs, &export_partitions);
861
862 // Zoltan_Destroy(&zz);
863
864 // liberation des tableau temporaires
865 freeConstraints();
866
867 cells_new_owner.synchronize();
869}
870
871/*---------------------------------------------------------------------------*/
872/*---------------------------------------------------------------------------*/
873
878
879/*---------------------------------------------------------------------------*/
880/*---------------------------------------------------------------------------*/
881
886ARCANE_REGISTER_SERVICE_ZOLTANMESHPARTITIONER(Zoltan,ZoltanMeshPartitioner);
887
888#if ARCANE_DEFAULT_PARTITIONER == ZOLTAN_DEFAULT_PARTITIONER
890 ServiceProperty("DefaultPartitioner",ST_SubDomain),
893ARCANE_REGISTER_SERVICE_ZOLTANMESHPARTITIONER(DefaultPartitioner,ZoltanMeshPartitioner);
894#endif
895
896/*---------------------------------------------------------------------------*/
897/*---------------------------------------------------------------------------*/
898
899} // End namespace Arcane
900
901/*---------------------------------------------------------------------------*/
902/*---------------------------------------------------------------------------*/
Types et macros pour itérer sur les entités du maillage.
#define ENUMERATE_FACE(name, group)
Enumérateur générique d'un groupe de faces.
#define ENUMERATE_CELL(name, group)
Enumérateur générique d'un groupe de mailles.
#define ENUMERATE_ITEM(name, group)
Enumérateur générique d'un groupe de noeuds.
#define ARCANE_SERVICE_INTERFACE(ainterface)
Macro pour déclarer une interface lors de l'enregistrement d'un service.
Integer size() const
Nombre d'éléments du vecteur.
CaseOptionsZoltanMeshPartitioner * options() const
Options du jeu de données du service.
ArcaneZoltanMeshPartitionerObject(const Arcane::ServiceBuildInfo &sbi)
Constructeur.
Exception lorsqu'un argument est invalide.
Vue modifiable d'un tableau d'un type T.
void copy(const U &copy_array)
Recopie le tableau copy_array dans l'instance.
void fill(const DataType &data)
Remplissage du tableau.
void resize(Int64 s)
Change le nombre d'éléments du tableau à s.
void reserve(Int64 new_capacity)
Réserve le mémoire pour new_capacity éléments.
Maille d'un maillage.
Definition Item.h:1214
Face d'une maille.
Definition Item.h:964
Cell frontCell() const
Maille devant la face (maille nulle si aucune)
Definition Item.h:1652
Cell backCell() const
Maille derrière la face (maille nulle si aucune)
Definition Item.h:1646
virtual ITraceMng * traceMng()=0
Gestionnaire de message associé
Interface d'un partitionneur de maillage.
Interface d'un partitionneur de maillage.
virtual IMesh * mesh() const =0
Maillage associé au partitionneur.
virtual VariableNodeReal3 & nodesCoordinates()=0
Coordonnées des noeuds.
Interface du gestionnaire de parallélisme pour un sous-domaine.
virtual Int32 commRank() const =0
Rang de cette instance dans le communicateur.
virtual void scan(eReduceType rt, ArrayView< char > v)=0
Applique un algorithme de prefix-um sur les valeurs de v via l'opération rt.
virtual Int32 commSize() const =0
Nombre d'instance dans le communicateur.
virtual void barrier()=0
Effectue une barière.
ItemUniqueId uniqueId() const
Identifiant unique sur tous les domaines.
Definition Item.h:225
constexpr bool null() const
true si l'entité est nul (i.e. non connecté au maillage)
Definition Item.h:216
Classe de base d'un service d'équilibrage de charge.
IMesh * mesh() const override
Maillage associé au partitionneur.
virtual void changeOwnersFromCells()
Positionne les nouveaux propriétaires des noeuds, arêtes et faces à partir des mailles.
Variable scalaire sur un type d'entité du maillage.
Noeud d'un maillage.
Definition Item.h:582
Classe gérant un vecteur de réel de dimension 3.
Definition Real3.h:132
Structure contenant les informations pour créer un service.
Propriétés de création d'un service.
Vecteur 1D de données avec sémantique par référence.
Chaîne de caractères unicode.
bool null() const
Retourne true si la chaîne est nulle.
Definition String.cc:305
TraceAccessor(ITraceMng *m)
Construit un accesseur via le gestionnaire de trace m.
TraceMessage fatal() const
Flot pour un message d'erreur fatale.
TraceMessage info() const
Flot pour un message d'information.
TraceMessage warning() const
Flot pour un message d'avertissement.
Vecteur 1D de données avec sémantique par valeur (style STL).
Paramètres nécessaires à la construction d'une variable.
Partitioneur de maillage utilisant la bibliotheque Zoltan.
virtual void notifyEndPartition()
Notification lors de la fin d'un re-partitionnement (après échange des entités)
virtual void build()
Construction de niveau build du service.
virtual void partitionMesh(bool initial_partition)
ItemGroupT< Cell > CellGroup
Groupe de mailles.
Definition ItemTypes.h:183
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro pour enregistrer un service.
MeshVariableScalarRefT< Node, Real3 > VariableNodeReal3
Grandeur au noeud de type coordonnées.
ItemVariableScalarRefT< Int32 > VariableItemInt32
Grandeur de type entier 32 bits.
Integer toInteger(Real r)
Converti un Real en Integer.
double toDouble(Real r)
Converti un Real en double.
@ ReduceSum
Somme des valeurs.
ARCCORE_BASE_EXPORT void sleep(Integer nb_second)
Met le process en sommeil pendant nb_second secondes.
ARCCORE_BASE_EXPORT String getEnvironmentVariable(const String &name)
Variable d'environnement du nom name.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
UniqueArray< Int64 > Int64UniqueArray
Tableau dynamique à une dimension d'entiers 64 bits.
Definition UtilsTypes.h:355
Int32 Integer
Type représentant un entier.
@ ST_SubDomain
Le service s'utilise au niveau du sous-domaine.
UniqueArray< Int32 > Int32UniqueArray
Tableau dynamique à une dimension d'entiers 32 bits.
Definition UtilsTypes.h:357
@ IK_Cell
Entité de maillage de genre maille.
double Real
Type représentant un réel.
std::int32_t Int32
Type entier signé sur 32 bits.