Arcane  v3.16.2.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
MetisMeshPartitioner.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/* MetisMeshPartitioner.cc (C) 2000-2025 */
9/* */
10/* Partitioneur de maillage utilisant la bibliothèque PARMetis. */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arcane/utils/StringBuilder.h"
15#include "arcane/utils/Iostream.h"
16#include "arcane/utils/PlatformUtils.h"
17#include "arcane/utils/Convert.h"
18#include "arcane/utils/NotImplementedException.h"
19#include "arcane/utils/ArgumentException.h"
20#include "arcane/utils/FloatingPointExceptionSentry.h"
21
22#include "arcane/ISubDomain.h"
23#include "arcane/IParallelMng.h"
24#include "arcane/ItemEnumerator.h"
25#include "arcane/IMesh.h"
26#include "arcane/IMeshSubMeshTransition.h"
27#include "arcane/ItemGroup.h"
28#include "arcane/Service.h"
29#include "arcane/Timer.h"
30#include "arcane/FactoryService.h"
31#include "arcane/ItemPrinter.h"
32#include "arcane/IItemFamily.h"
33#include "arcane/MeshVariable.h"
34#include "arcane/VariableBuildInfo.h"
35
36#include "arcane/std/MeshPartitionerBase.h"
37#include "arcane/std/MetisMeshPartitioner_axl.h"
38
39// Au cas où on utilise mpich2 ou openmpi pour éviter d'inclure le C++
40#define MPICH_SKIP_MPICXX
41#define OMPI_SKIP_MPICXX
42#include <parmetis.h>
43// #define CEDRIC_LB_2
44
45#if (PARMETIS_MAJOR_VERSION < 4)
46#error "Your version of Parmetis is too old .Version 4.0+ is required"
47#endif
48
49typedef real_t realtype;
50typedef idx_t idxtype;
51
52#include "arcane/std/PartitionConverter.h"
53#include "arcane/std/GraphDistributor.h"
54#include "arcane/std/internal/MetisWrapper.h"
55
56#include <chrono>
57
58/*---------------------------------------------------------------------------*/
59/*---------------------------------------------------------------------------*/
60
61namespace Arcane
62{
63
64using MetisCallStrategy = TypesMetisMeshPartitioner::MetisCallStrategy;
65using MetisEmptyPartitionStrategy = TypesMetisMeshPartitioner::MetisEmptyPartitionStrategy;
66
67/*---------------------------------------------------------------------------*/
68/*---------------------------------------------------------------------------*/
72class MetisMeshPartitioner
74{
75 public:
76
77 explicit MetisMeshPartitioner(const ServiceBuildInfo& sbi);
78
79 public:
80
81 void build() override {}
82
83 public:
84
85 void partitionMesh(bool initial_partition) override;
86 void partitionMesh(bool initial_partition,Int32 nb_part) override;
87
88
89 private:
90
91 IParallelMng* m_parallel_mng = nullptr;
92 Integer m_nb_refine = -1;
93 Integer m_random_seed = 15;
94 bool m_disable_floatingexception = false;
95 void _partitionMesh(bool initial_partition,Int32 nb_part);
96 void _removeEmptyPartsV1(Int32 nb_part, Int32 nb_own_cell, ArrayView<idxtype> metis_part);
97 void _removeEmptyPartsV2(Int32 nb_part,ArrayView<idxtype> metis_part);
98 Int32 _removeEmptyPartsV2Helper(Int32 nb_part,ArrayView<idxtype> metis_part,Int32 algo_iteration);
100 Span<const idxtype> metis_vtkdist,
101 Span<const idxtype> metis_xadj,
102 Span<const idxtype> metis_adjncy,
103 Span<const idxtype> metis_vwgt,
104 Span<const idxtype> metis_ewgt,
105 const String& Name) const;
106};
107
108/*---------------------------------------------------------------------------*/
109/*---------------------------------------------------------------------------*/
110
111MetisMeshPartitioner::
112MetisMeshPartitioner(const ServiceBuildInfo& sbi)
113: ArcaneMetisMeshPartitionerObject(sbi)
114{
115 m_parallel_mng = mesh()->parallelMng();
116 String s = platform::getEnvironmentVariable("ARCANE_DISABLE_METIS_FPE");
117 if (s=="1" || s=="TRUE")
118 m_disable_floatingexception = true;
119}
120
121/*---------------------------------------------------------------------------*/
122/*---------------------------------------------------------------------------*/
123
125partitionMesh(bool initial_partition)
126{
127 Int32 nb_part = m_parallel_mng->commSize();
128 partitionMesh(initial_partition,nb_part);
129}
130
131/*---------------------------------------------------------------------------*/
132/*---------------------------------------------------------------------------*/
133
135partitionMesh(bool initial_partition,Int32 nb_part)
136{
137 // Signale que le partitionnement peut planter, car metis n'est pas toujours
138 // tres fiable sur le calcul flottant.
139 initial_partition = (subDomain()->commonVariables().globalIteration() <= 2);
140 if (m_disable_floatingexception){
142 _partitionMesh(initial_partition,nb_part);
143 }
144 else
145 _partitionMesh(initial_partition,nb_part);
146}
147
148/*---------------------------------------------------------------------------*/
149/*---------------------------------------------------------------------------*/
150
152_partitionMesh(bool initial_partition,Int32 nb_part)
153{
154 ISubDomain* sd = subDomain();
155 IParallelMng* pm = m_parallel_mng;
156 Int32 my_rank = pm->commRank();
157 Int32 nb_rank = pm->commSize();
158 IMesh* mesh = this->mesh();
159 bool force_partition = false;
160
161 if (nb_part<nb_rank)
162 ARCANE_THROW(ArgumentException,"partition with nb_part ({0}) < nb_rank ({1})",
163 nb_part,nb_rank);
164
165 if (nb_rank==1){
166 info() << "INFO: Unable to test load balancing on a single sub-domain";
167 return;
168 }
169
170 // TODO : comprendre les modifs de l'IFPEN
171 // utilise toujours re-partitionnement complet.
172 // info() << "Metis: params " << m_nb_refine << " " << imbalance() << " " << maxImbalance();
173 // info() << "Metis: params " << (m_nb_refine>=10) << " " << (imbalance()>4.0*maxImbalance()) << " " << (imbalance()>1.0);
174 bool force_full_repartition = false;
175 //force_full_repartition = true;
176 force_full_repartition |= (m_nb_refine < 0); // toujours la première fois, pour compatibilité avec l'ancienne version
177
178 if (nb_part != nb_rank) { // Pas de "repartitionnement" si le nombre de parties ne correspond pas au nombre de processeur.
179 force_full_repartition = true;
180 force_partition = true;
181 }
182
183 info() << "WARNING: compensating the potential lack of 'Metis' options in case of manual instanciation";
184 Integer max_diffusion_count = 10; // reprise des options par défaut du axl
185 Real imbalance_relative_tolerance = 4;
186 float tolerance_target = 1.05f;
187 bool dump_graph = false;
188 MetisCallStrategy call_strategy = MetisCallStrategy::one_processor_per_node;
189 bool in_out_digest = false;
190
191 // Variables d'environnement permettant de regler les options lorsqu'il n'y a pas de jeu de donnees
192
193 String call_strategy_env = platform::getEnvironmentVariable("ARCANE_METIS_CALL_STRATEGY");
194 if (call_strategy_env == "all-processors"){
195 call_strategy = MetisCallStrategy::all_processors;
196 } else if (call_strategy_env == "one-processor-per-node") {
197 call_strategy = MetisCallStrategy::one_processor_per_node;
198 } else if (call_strategy_env == "two-processors-two-nodes") {
199 call_strategy = MetisCallStrategy::two_processors_two_nodes;
200 } else if (call_strategy_env == "two-gathered-processors") {
201 call_strategy = MetisCallStrategy::two_gathered_processors;
202 } else if (call_strategy_env == "two-scattered-processors") {
203 call_strategy = MetisCallStrategy::two_scattered_processors;
204 } else if (!call_strategy_env.null()) {
205 ARCANE_FATAL("Invalid value '{0}' for ARCANE_METIS_CALL_STRATEGY environment variable",call_strategy_env);
206 }
207
208 if (auto v = Convert::Type<Int32>::tryParseFromEnvironment("ARCANE_METIS_INPUT_OUTPUT_DIGEST", true))
209 in_out_digest = (v.value()!=0);
210 if (auto v = Convert::Type<Int32>::tryParseFromEnvironment("ARCANE_METIS_DUMP_GRAPH", true))
211 dump_graph = (v.value()!=0);
212
213 if (options()){
214 max_diffusion_count = options()->maxDiffusiveCount();
215 imbalance_relative_tolerance = options()->imbalanceRelativeTolerance();
216 tolerance_target = (float)(options()->toleranceTarget());
217 dump_graph = options()->dumpGraph();
218 call_strategy = options()->metisCallStrategy();
219 in_out_digest = options()->inputOutputDigest();
220 }
221
222 if (max_diffusion_count > 0)
223 force_full_repartition |= (m_nb_refine>=max_diffusion_count);
224 force_full_repartition |= (imbalance()>imbalance_relative_tolerance*maxImbalance());
225 force_full_repartition |= (imbalance()>1.0);
226
227 // initialisations pour la gestion des contraintes (sauf initUidRef)
228 initConstraints(false);
229
230 bool is_shared_memory = pm->isThreadImplementation();
231 // En mode mémoire partagé, pour l'instant on force le fait d'utiliser un seul
232 // processeur par noeud car on ne sait pas si on dispose de plusieurs rang par noeud.
233 // Notamment, en mode mémoire partagé sans MPI on n'a obligatoirement qu'un seul noeud
234 if (is_shared_memory)
235 call_strategy = MetisCallStrategy::one_processor_per_node;
236
237 idxtype nb_weight = nbCellWeight();
238
239 info() << "Load balancing with Metis nb_weight=" << nb_weight << " initial=" << initial_partition
240 << " call_strategy=" << (int)call_strategy
241 << " is_shared_memory?=" << is_shared_memory
242 << " disabling_fpe?=" << m_disable_floatingexception
243 << " sizeof(idxtype)==" << sizeof(idxtype);
244
245 if (nb_weight==0)
246 initial_partition = true;
247
248 UniqueArray<idxtype> metis_vtkdist(nb_rank+1);
249 Integer total_nb_cell = 0;
250
251 // Contient les numéros uniques des entités dans la renumérotation
252 // propre à metis
253 VariableCellInteger cell_metis_uid(VariableBuildInfo(mesh,"CellsMetisUid",IVariable::PNoDump));
254
255 UniqueArray<Integer> global_nb_own_cell(nb_rank);
256 CellGroup own_cells = mesh->ownCells();
257 Integer nb_own_cell = nbOwnCellsWithConstraints(); // on tient compte des contraintes
258 pm->allGather(ConstArrayView<Integer>(1,&nb_own_cell),global_nb_own_cell);
259 Int32 nb_empty_part = 0;
260 {
261 //Integer total_first_uid = 0;
262 metis_vtkdist[0] = 0;
263 for( Integer i=0; i<nb_rank; ++i ){
264 //total_first_uid += global_nb_own_cell[i];
265 Int32 n = global_nb_own_cell[i];
266 if (n==0)
267 ++nb_empty_part;
268 total_nb_cell += n;
269 metis_vtkdist[i+1] = static_cast<idxtype>(total_nb_cell);
270 // info() << "METIS VTKDIST " << (i+1) << ' ' << metis_vtkdist[i+1];
271 }
272 }
273 // N'appelle pas Parmetis si on n'a pas de mailles car sinon cela provoque une erreur.
274 info() << "Total nb_cell=" << total_nb_cell << " nb_empty_partition=" << nb_empty_part;
275 if (total_nb_cell==0){
276 info() << "INFO: mesh '" << mesh->name() << " has no cell. No partitioning is needed";
277 freeConstraints();
278 return;
279 }
280
281 // HP: Skip this shortcut because it produces undefined itemsNewOwner variable.
282/*
283 if (total_nb_cell < nb_rank) { // Dans ce cas, pas la peine d'appeler ParMetis.
284 warning() << "There are no subdomain except cells, no reffinement";
285 freeConstraints();
286 return;
287 }
288*/
289
290 // Nombre max de mailles voisines connectées aux mailles
291 // en supposant les mailles connectées uniquement par les faces
292 // (ou les arêtes pour une maille 2D dans un maillage 3D)
293 // Cette valeur sert à préallouer la mémoire pour la liste des mailles voisines
294 Integer nb_max_face_neighbour_cell = 0;
295 {
296 // Renumérote les mailles pour Metis pour que chaque sous-domaine
297 // ait des mailles de numéros consécutifs
298 Integer mid = static_cast<Integer>(metis_vtkdist[my_rank]);
299 ENUMERATE_ (Cell, i_item, own_cells) {
300 Cell item = *i_item;
301 if (cellUsedWithConstraints(item)){
302 cell_metis_uid[item] = mid;
303 ++mid;
304 }
305 if (item.hasFlags(ItemFlags::II_HasEdgeFor1DItems))
306 nb_max_face_neighbour_cell += item.nbEdge();
307 else
308 nb_max_face_neighbour_cell += item.nbFace();
309 }
310 cell_metis_uid.synchronize();
311 }
312
313 _initUidRef(cell_metis_uid);
314
315 // libération mémoire
316 cell_metis_uid.setUsed(false);
317
318 SharedArray<idxtype> metis_xadj;
319 metis_xadj.reserve(nb_own_cell+1);
320
321 SharedArray<idxtype> metis_adjncy;
322 metis_adjncy.reserve(nb_max_face_neighbour_cell);
323
324
325 // Construction de la connectivité entre les cellules et leurs voisines en tenant compte des contraintes
326 // (La connectivité se fait suivant les faces)
327 Int64UniqueArray neighbour_cells;
328 UniqueArray<float> edge_weights;
329 edge_weights.resize(0);
330 ENUMERATE_CELL(i_item,own_cells){
331 Cell item = *i_item;
332
333 if (!cellUsedWithConstraints(item))
334 continue;
335
336 metis_xadj.add(metis_adjncy.size());
337
338 getNeighbourCellsUidWithConstraints(item, neighbour_cells, &edge_weights);
339 for( Integer z=0; z<neighbour_cells.size(); ++z )
340 metis_adjncy.add(static_cast<idxtype>(neighbour_cells[z]));
341 }
342 metis_xadj.add(metis_adjncy.size());
343
344 idxtype wgtflag = 3; // Sommets et aretes
345//2; // Indique uniquement des poids sur les sommets
346 idxtype numflags = 0;
347 idxtype nparts = static_cast<int>(nb_part);
348
349 // Poids aux sommets du graphe (les mailles)
350// if (initial_partition)
351// nb_weight = 1;
352
353 String s = platform::getEnvironmentVariable("ARCANE_LB_PARAM");
354 if (!s.null() && (s.upper() == "PARTICLES"))
355 nb_weight = 1; // Only one weight in this case !
356
357#if CEDRIC_LB_2
358 nb_weight = 1;
359#endif
360 bool cond = (nb_weight == 1);
361 // Real max_int = 1 << (sizeof(idxtype)*8-1);
362 UniqueArray<realtype> metis_ubvec(CheckedConvert::toInteger(nb_weight));
363 SharedArray<float> cells_weights;
364
365 PartitionConverter<float,idxtype> converter(pm, (Real)IDX_MAX, cond);
366
367
373 cells_weights = cellsWeightsWithConstraints(CheckedConvert::toInteger(nb_weight));
374 // Déséquilibre autorisé pour chaque contrainte
375 metis_ubvec.resize(CheckedConvert::toInteger(nb_weight));
376 metis_ubvec.fill(tolerance_target);
377
378 converter.reset(CheckedConvert::toInteger(nb_weight));
379
380 do {
381 cond = converter.isBalancable<realtype>(cells_weights, metis_ubvec, nb_part);
382 if (!cond) {
383 info() << "We have to increase imbalance :";
384 info() << metis_ubvec;
385 for (auto& tol: metis_ubvec ) {
386 tol *= 1.5f;
387 }
388 }
389 } while (!cond);
390
391 ArrayConverter<float,idxtype,PartitionConverter<float,idxtype> > metis_vwgtConvert(cells_weights, converter);
392
393 SharedArray<idxtype> metis_vwgt(metis_vwgtConvert.array().constView());
394
395 const bool do_print_weight = false;
396 if (do_print_weight){
397 StringBuilder fname;
398 Integer iteration = mesh->subDomain()->commonVariables().globalIteration();
399 fname = "weigth-";
400 fname += pm->commRank();
401 fname += "-";
402 fname += iteration;
403 String f(fname);
404 std::ofstream dumpy(f.localstr());
405 for (int i = 0; i < metis_xadj.size()-1 ; ++i) {
406 dumpy << " Weight uid=" << i;
407 for( int j=0 ; j < nb_weight ; ++j ){
408 dumpy << " w=" << *(metis_vwgt.begin()+i*nb_weight+j)
409 << "(" << cells_weights[i*nb_weight+j] <<")";
410 }
411 dumpy << '\n';
412 }
413 }
414
415 converter.reset(1); // Only one weight for communications !
416
417 converter.computeContrib(edge_weights);
418 ArrayConverter<float,idxtype,PartitionConverter<float,idxtype> > metis_ewgt_convert(edge_weights, converter);
419 SharedArray<idxtype> metis_ewgt((UniqueArray<idxtype>)metis_ewgt_convert.array());
420
421 SharedArray<float> cells_size;
422 realtype itr=50; // In average lb every 20 iterations ?
423 SharedArray<idxtype> metis_vsize;
424 if (!initial_partition) {
425 cells_size = cellsSizeWithConstraints();
426 converter.computeContrib(cells_size, (Real)itr);
427 metis_vsize.resize(metis_xadj.size()-1,1);
428 }
429
430
431 idxtype metis_options[4];
432 metis_options[0] = 0;
433
434 // By default RandomSeed is fixed
435
436#ifdef CEDRIC_LB_2
437 m_random_seed = Convert::toInteger((Real)MPI_Wtime());
438#endif // CEDRIC_LB_2
439
440 metis_options[0] = 1;
441 metis_options[1] = 0;
442 metis_options[2] = m_random_seed;
443 metis_options[3] = 1; // Initial distribution information is implicit
444 /*
445 * Comment to keep same initial random seed each time !
446 m_random_seed++;
447 */
448
449 idxtype metis_edgecut = 0;
450
451 // TODO: compute correct part number !
452 UniqueArray<idxtype> metis_part;
453
454 std::chrono::high_resolution_clock clock;
455 auto start_time = clock.now();
456
457 GraphDistributor gd(pm);
458
459 // Il y a actuellement 2 mecanismes de regroupement du graph : celui fait par "GraphDistributor" et
460 // celui fait par le wrapper ParMetis. Il est possible de combiner les 2 (two_processors_two_nodes).
461 // A terme, ces 2 mecanismes devraient fusionner et il ne faudrait conserver que le "GraphDistributor".
462 //
463 // La redistribution par le wrapper n'est faite que dans les 2 cas suivants :
464 // - on veut un regroupement sur 2 processeurs apres un regroupement par noeuds (two_processors_two_nodes)
465 // - on veut un regroupement direct sur 2 processeurs, en esperant qu'ils soient sur 2 noeuds distincts (two_scattered_processors)
466
467 bool redistribute = true; // redistribution par GraphDistributor
468 bool redistribute_by_wrapper = false; // redistribution par le wrapper ParMetis
469
470 if (call_strategy == MetisCallStrategy::all_processors || call_strategy == MetisCallStrategy::two_scattered_processors) {
471 redistribute = false;
472 }
473
474 if (call_strategy == MetisCallStrategy::two_processors_two_nodes || call_strategy == MetisCallStrategy::two_scattered_processors) {
475 redistribute_by_wrapper = true;
476 }
477 // Indique si on n'autorise de n'utiliser qu'un seul PE.
478 // historiquement (avant avril 2020) cela n'était pas autorisé car cela revenait
479 // à appeler 'ParMetis' sur un seul PE ce qui n'était pas supporté.
480 // Maintenant, on appelle directement Metis dans ce cas donc cela ne pose pas de
481 // problèmes. Cependant pour des raisons historiques on garde l'ancien comportement
482 // sauf pour deux cas:
483 // - si l'échange de messages est en mode mémoire partagé. Comme le partitionnement
484 // dans ce mode n'était pas supporté avant, il n'y a
485 // pas d'historique à conserver. De plus cela est indispensable si on n'utilise
486 // qu'un seul noeud de calcul car alors
487 // - lors du partionnement initial et s'il y a des partitions vides. Ce cas n'existait
488 // pas avant non plus car on utilisait 'MeshPartitionerTester' pour faire un
489 // premier pré-partitionnement qui garantit aucune partition vide. Cela permet
490 // d'éviter ce pré-partitionnement.
491
492 bool gd_allow_only_one_rank = false;
493 if (is_shared_memory || (nb_empty_part!=0 && initial_partition))
494 gd_allow_only_one_rank = true;
495
496 // S'il n'y a qu'une seule partie non-vide on n'utilise qu'un seul rang. Ce cas
497 // intervient normalement uniquement en cas de partitionnement initial et si
498 // un seul rang (en général le rang 0) a des mailles.
499 if ((nb_empty_part+1)==nb_rank){
500 info() << "Initialize GraphDistributor with max rank=1";
501 gd.initWithMaxRank(1);
502 }
503 else if (call_strategy == MetisCallStrategy::two_gathered_processors && (nb_rank > 2)) {
504 // Seuls les 2 premiers processeurs seront utilisés.
505 info() << "Initialize GraphDistributor with max rank=2";
506 gd.initWithMaxRank(2);
507 }
508 else {
509 info() << "Initialize GraphDistributor with one rank per node"
510 << " (allow_one?=" << gd_allow_only_one_rank << ")";
511 gd.initWithOneRankPerNode(gd_allow_only_one_rank);
512 }
513
514 IParallelMng* metis_pm = pm;
515
516 info() << "Using GraphDistributor?=" << redistribute;
517 if (redistribute){
518 metis_xadj = gd.convert<idxtype>(metis_xadj, &metis_part, true);
519 metis_vwgt = gd.convert<idxtype>(metis_vwgt);
520 metis_adjncy = gd.convert<idxtype>(metis_adjncy);
521 metis_ewgt = gd.convert<idxtype>(metis_ewgt);
522 if (!initial_partition)
523 metis_vsize = gd.convert<idxtype>(metis_vsize);
524 metis_pm = gd.subParallelMng();
525
526 metis_vtkdist.resize(gd.size()+1);
527 if (gd.contribute()) {
528 UniqueArray<Integer> buffer(gd.size());
529 Integer nbVertices = metis_part.size();
530 gd.subParallelMng()->allGather(ConstArrayView<Integer>(1,&nbVertices),buffer);
531 metis_vtkdist[0] = 0;
532 for (int i = 1 ; i < metis_vtkdist.size() ; ++i) {
533 metis_vtkdist[i] = metis_vtkdist[i-1] + buffer[i-1];
534 }
535 }
536 metis_options[3] = PARMETIS_PSR_UNCOUPLED ; // Initial distribution information is in metis_part
537 }
538 else{
539
540 metis_part = UniqueArray<idxtype>(nb_own_cell, my_rank);
541 }
542 MPI_Comm metis_mpicomm = static_cast<MPI_Comm>(metis_pm->communicator());
543
544 bool do_call_metis = true;
545 if (redistribute)
546 do_call_metis = gd.contribute();
547
548 if (do_call_metis) {
549 if (dump_graph) {
550 Integer iteration = sd->commonVariables().globalIteration();
551 StringBuilder name("graph-");
552 name += iteration;
553 _writeGraph(metis_pm, metis_vtkdist, metis_xadj, metis_adjncy, metis_vwgt, metis_ewgt, name.toString());
554 }
555
556 // ParMetis >= 4 requires to define tpwgt.
557 const Integer tpwgt_size = CheckedConvert::toInteger(nb_part) * CheckedConvert::toInteger(nb_weight);
558 UniqueArray<realtype> tpwgt(tpwgt_size);
559 float fill_value = (float)(1.0/(double)nparts);
560 tpwgt.fill(fill_value);
561
562 int retval = METIS_OK;
563 Timer::Action ts(sd,"Metis");
564
565 MetisWrapper wrapper(metis_pm);
566 force_partition |= initial_partition;
567 force_full_repartition |= force_partition;
568 if (force_full_repartition){
569 m_nb_refine = 0;
570 if (force_partition) {
571 info() << "Metis: use a complete partitioning.";
572 retval = wrapper.callPartKway(in_out_digest,
573 redistribute_by_wrapper,
574 metis_vtkdist.data(),
575 metis_xadj.data(),
576 metis_adjncy.data(),
577 metis_vwgt.data(),
578 metis_ewgt.data(),
579 &wgtflag,&numflags,&nb_weight,
580 &nparts, tpwgt.data(),
581 metis_ubvec.data(),
582 metis_options,&metis_edgecut,
583 metis_part.data());
584 }
585 else {
586 info() << "Metis: use a complete REpartitioning.";
587 retval = wrapper.callAdaptiveRepart(in_out_digest,
588 redistribute_by_wrapper,
589 metis_vtkdist.data(),
590 metis_xadj.data(),
591 metis_adjncy.data(),
592 metis_vwgt.data(),
593 metis_vsize.data(), // Vsize !
594 metis_ewgt.data(),
595 &wgtflag,&numflags,&nb_weight,
596 &nparts,tpwgt.data(),
597 metis_ubvec.data(),
598 &itr, metis_options,&metis_edgecut,
599 metis_part.data());
600 }
601 }
602 else{
603 // TODO: mettre dans MetisWrapper et supprimer utilisation de .data()
604 // (cela doit être faire par le wrapper)
605 ++m_nb_refine;
606 info() << "Metis: use a diffusive REpartitioning";
607 retval = ParMETIS_V3_RefineKway(metis_vtkdist.data(),
608 metis_xadj.data(),
609 metis_adjncy.data(),
610 metis_vwgt.data(),
611 metis_ewgt.data(),
612 &wgtflag,&numflags,&nb_weight,
613 &nparts,tpwgt.data(),
614 metis_ubvec.data(),
615 metis_options,&metis_edgecut,
616 metis_part.data(),
617 &metis_mpicomm);
618 }
619 if (retval != METIS_OK) {
620 ARCANE_FATAL("Error while computing ParMetis partitioning r='{0}'",retval);
621 }
622 }
623 if (redistribute){
624 metis_part = gd.convertBack<idxtype>(metis_part, nb_own_cell);
625 }
626
627 auto end_time = clock.now();
628 std::chrono::microseconds duration = std::chrono::duration_cast<std::chrono::microseconds>(end_time - start_time);
629 info() << "Time to partition using parmetis = " << duration.count() << " us ";
630
631 // Stratégie à adopter pour supprimer les partitions vides.
632 // A noter qu'il faut faire cela avant d'appliquer les éventuelles contraintes
633 // pour garantir que ces dernières seront bien respectées
634 if (options()){
635 switch(options()->emptyPartitionStrategy()){
636 case MetisEmptyPartitionStrategy::DoNothing:
637 break;
638 case MetisEmptyPartitionStrategy::TakeFromBiggestPartitionV1:
639 _removeEmptyPartsV1(nb_part,nb_own_cell,metis_part);
640 break;
641 case MetisEmptyPartitionStrategy::TakeFromBiggestPartitionV2:
642 _removeEmptyPartsV2(nb_part,metis_part);
643 break;
644 }
645 }
646 else
647 _removeEmptyPartsV2(nb_part,metis_part);
648
649 VariableItemInt32& cells_new_owner = mesh->toPrimaryMesh()->itemsNewOwner(IK_Cell);
650 {
651 Integer index = 0;
652 ENUMERATE_CELL(i_item,own_cells){
653 Cell item = *i_item;
654 if (!cellUsedWithConstraints(item))
655 continue;
656
657 Int32 new_owner = CheckedConvert::toInt32(metis_part[index]);
658 ++index;
659 changeCellOwner(item, cells_new_owner, new_owner);
660 }
661 }
662
663 // libération des tableaux temporaires
664 freeConstraints();
665
666 cells_new_owner.synchronize();
667
669
670 m_nb_refine = m_parallel_mng->reduce(Parallel::ReduceMax, m_nb_refine);
671}
672
673/*---------------------------------------------------------------------------*/
674/*---------------------------------------------------------------------------*/
686_removeEmptyPartsV1(const Int32 nb_part,const Int32 nb_own_cell,ArrayView<idxtype> metis_part)
687{
688 IParallelMng* pm = m_parallel_mng;
689 Int32 my_rank = pm->commRank();
690
691 //The following code insures that no empty part will be created.
692 // TODO: faire une meilleure solution, qui prend des elements sur la partie la plus lourde.
693 UniqueArray<Int32> elem_by_part(nb_part,0);
694 UniqueArray<Int32> min_part(nb_part);
695 UniqueArray<Int32> max_part(nb_part);
696 UniqueArray<Int32> sum_part(nb_part);
697 UniqueArray<Int32> min_roc(nb_part);
698 UniqueArray<Int32> max_proc(nb_part);
699 for (int i =0 ; i < nb_own_cell ; i++) {
700 elem_by_part[CheckedConvert::toInteger(metis_part[i])]++;
701 }
702 pm->computeMinMaxSum(elem_by_part, min_part, max_part, sum_part, min_roc, max_proc);
703
704 int nb_hole=0;
705 Int32 max_part_id = -1;
706 Int32 max_part_nbr = -1;
707 // Compute number of empty parts
708 for(int i = 0; i < nb_part ; i++) {
709 if (sum_part[i] == 0) {
710 nb_hole ++;
711 }
712 if(max_part_nbr < sum_part[i]) {
713 max_part_nbr = sum_part[i];
714 max_part_id = i;
715 }
716 }
717 info() << "Parmetis: number empty parts " << nb_hole;
718
719 // Le processeur ayant la plus grosse partie de la plus
720 // grosse partition comble les trous
721 if(my_rank == max_proc[max_part_id]) {
722 int offset = 0;
723 for(int i = 0; i < nb_part ; i++) {
724 // On ne comble que s'il reste des mailles
725 if (sum_part[i] == 0 && offset < nb_own_cell) {
726 while(offset < nb_own_cell && metis_part[offset] != max_part_id) {
727 offset++;
728 }
729 // Si on est sorti du while car pas assez de mailles
730 if(offset == nb_own_cell)
731 break;
732 // Le trou est comblé par ajout d'une seule maille
733 metis_part[offset] = i;
734 offset++;
735 }
736 }
737 }
738}
739
740/*---------------------------------------------------------------------------*/
741/*---------------------------------------------------------------------------*/
750_removeEmptyPartsV2Helper(const Int32 nb_part,ArrayView<idxtype> metis_part,Int32 algo_iteration)
751{
752 IParallelMng* pm = m_parallel_mng;
753 Int32 my_rank = pm->commRank();
754 const Int32 nb_own_cell = metis_part.size();
755 // The following code insures that no empty part will be created.
756 UniqueArray<idxtype> elem_by_part(nb_part,0);
757 UniqueArray<idxtype> min_part(nb_part);
758 UniqueArray<idxtype> max_part(nb_part);
759 UniqueArray<idxtype> sum_part(nb_part);
760 UniqueArray<Int32> min_rank(nb_part);
761 UniqueArray<Int32> max_rank(nb_part);
762 for (int i =0 ; i < nb_own_cell ; i++) {
763 elem_by_part[metis_part[i]]++;
764 }
765 pm->computeMinMaxSum(elem_by_part, min_part, max_part, sum_part, min_rank, max_rank);
766
767 // Rang des parties vides
768 UniqueArray<Int32> empty_part_ranks;
769 int nb_hole = 0;
770 Int32 max_part_id = -1;
771 Int64 max_part_nbr = -1;
772 // Compute number of empty parts
773 Int64 total_nb_cell = 0;
774 for(int i = 0; i < nb_part ; i++) {
775 Int64 current_sum_part = sum_part[i];
776 total_nb_cell += current_sum_part;
777 if (current_sum_part == 0) {
778 empty_part_ranks.add(i);
779 nb_hole ++;
780 }
781 if (max_part_nbr < current_sum_part) {
782 max_part_nbr = current_sum_part;
783 max_part_id = i;
784 }
785 }
786 if (max_part_id<0)
787 ARCANE_FATAL("Bad value max_part_id ({0})",max_part_id);
788 info() << "Parmetis: check empty parts: (" << algo_iteration << ") nb_empty_parts=" << nb_hole
789 << " nb_max_part=" << max_part_nbr
790 << " max_part_rank=" << max_part_id
791 << " max_proc_max_part_id=" << max_rank[max_part_id]
792 << " empty_part_ranks=" << empty_part_ranks
793 << " total_nb_cell=" << total_nb_cell;
794
795 if (nb_hole==0)
796 return 0;
797
798 // Le processeur ayant la plus grosse partie de la plus
799 // grosse partition comble les trous
800 if (my_rank == max_rank[max_part_id]) {
801 // On garantit qu'on n'enlèvera pas toutes nos mailles
802 Int32 max_remove_cell = nb_own_cell - 1;
803 int offset = 0;
804 for(int i = 0; i < nb_part ; i++) {
805 // On ne comble que s'il reste des mailles
806 if (sum_part[i] == 0 && offset < nb_own_cell) {
807 while (offset < max_remove_cell && metis_part[offset] != max_part_id) {
808 offset++;
809 }
810 // Si on est sorti du while car pas assez de mailles
811 if (offset == max_remove_cell)
812 break;
813 // Le trou est comblé par ajout d'une seule maille
814 metis_part[offset] = i;
815 offset++;
816 }
817 }
818 }
819
820 return nb_hole;
821}
822
823/*---------------------------------------------------------------------------*/
824/*---------------------------------------------------------------------------*/
832_removeEmptyPartsV2(const Int32 nb_part,ArrayView<idxtype> metis_part)
833{
834 bool do_continue = true;
835 Int32 last_nb_hole = 0;
836 while (do_continue){
837 Int32 nb_hole = _removeEmptyPartsV2Helper(nb_part,metis_part,0);
838 info() << "Parmetis: nb_empty_partition=" << nb_hole << " last_nb_partition=" << last_nb_hole;
839 if (nb_hole==0)
840 break;
841 // Garanti qu'on sort de la boucle si on a toujours le même nombre de trous
842 // qu'avant. Cela permet d'éviter les boucles infinies.
843 // Cela signifie aussi qu'on ne peut pas combler les trous et donc il y
844 // aura surement des partitions vides
845 if (last_nb_hole>0 && last_nb_hole<=nb_hole){
846 pwarning() << "Can not remove all empty partitions. This is probably because you try"
847 << " to cut in too many partitions";
848 break;
849 }
850 last_nb_hole = nb_hole;
851 }
852}
853
854/*---------------------------------------------------------------------------*/
855/*---------------------------------------------------------------------------*/
856
862 Span<const idxtype> metis_vtkdist,
863 Span<const idxtype> metis_xadj,
864 Span<const idxtype> metis_adjncy,
865 Span<const idxtype> metis_vwgt,
866 Span<const idxtype> metis_ewgt,
867 const String& name) const
868{
869 int retval = 0;
870
871 idxtype nvtx = 0;
872 idxtype nwgt = 0;
873 bool have_ewgt = false;
874
875 Int32 commrank = pm->commRank();
876 Int32 commsize = pm->commSize();
877 info() << "COMM_SIZE=" << commsize << " RANK=" << commrank;
878 traceMng()->flush();
879 //MPI_Comm_size(metis_comm, &commsize);
880 //MPI_Comm_rank(metis_comm ,&commrank);
881 // NOTE GG: la gestion des erreurs via METIS_ERROR ne fonctionne pas: cela produit un blocage
882 // car il y a un MPI_Allreduce dans la partie sans erreur.
883 // NOTE GG: Ne pas utiliser MPI directement.
884
885 info() << "MetisVtkDist=" << metis_vtkdist;
886 info() << "MetisXAdj =" << metis_xadj;
887 info() << "MetisAdjncy =" << metis_adjncy;
888 info() << "MetisVWgt =" << metis_vwgt;
889 info() << "MetisEWgt =" << metis_ewgt;
890
891#define METIS_ERROR ARCANE_FATAL("_writeGraph")
892
893 StringBuilder filename(name);
894 filename += "_";
895 filename += commrank;
896 std::ofstream file(filename.toString().localstr());
897
898 if (metis_vtkdist.size() != commsize + 1)
899 METIS_ERROR;
900
901 nvtx = metis_vtkdist[commrank+1] - metis_vtkdist[commrank];
902 if (metis_xadj.size() != nvtx + 1) {
903 std::cerr << "_writeGraph : nvtx+1 = " << nvtx << " != " << metis_xadj.size() << std::endl;
904 METIS_ERROR;
905 }
906
907 if (nvtx != 0)
908 nwgt = metis_vwgt.size()/nvtx;
909 if (nwgt != 0 && metis_vwgt.size() % nwgt != 0)
910 METIS_ERROR;
911
912 have_ewgt = (metis_ewgt.size() != 0) ;
913 if (have_ewgt && metis_ewgt.size() != metis_adjncy.size())
914 METIS_ERROR;
915
916 if (!file.is_open())
917 METIS_ERROR;
918
919 Int64 localEdges = metis_xadj[nvtx];
920 Int64 globalEdges = pm->reduce(Parallel::ReduceSum,localEdges);
921
922 file << "2" << std::endl;
923 file << commsize << "\t" << commrank << std::endl;
924 file << metis_vtkdist[commsize] << "\t" << globalEdges << std::endl;
925 file << nvtx << "\t" << localEdges << std::endl;
926 file << "0\t" << "0" << have_ewgt << nwgt << std::endl;
927
928 /*
929 Each of these lines begins with the vertex label,
930 if necessary, the vertex load, if necessary, and the vertex degree, followed by the
931 description of the arcs. An arc is defined by the load of the edge, if necessary, and
932 by the label of its other end vertex.
933 */
934
935 //** Lbl Wgt Degree edWgt edLbl ...
936
937 for (idxtype vertnum = 0 ; vertnum < nvtx ; vertnum++) {
938 // file << vertnum + metis_vtkdist[commrank] << " ";
939 for (int dim = 0 ; dim < nwgt ; ++ dim)
940 file << metis_vwgt[vertnum*nwgt+dim] << '\t';
941 file << metis_xadj[vertnum + 1] - metis_xadj[vertnum];
942 for (idxtype edgenum = metis_xadj[vertnum] ;
943 edgenum < metis_xadj[vertnum + 1] ; ++edgenum) {
944 if (have_ewgt)
945 file << '\t' << metis_ewgt[edgenum];
946 file << '\t' << metis_adjncy[edgenum];
947 }
948 file << std::endl;
949 }
950 file.close();
951
952 pm->barrier();
953 return retval;
954}
955
956
957/*---------------------------------------------------------------------------*/
958/*---------------------------------------------------------------------------*/
959
964
965ARCANE_REGISTER_SERVICE_METISMESHPARTITIONER(Metis,MetisMeshPartitioner);
966
967#if ARCANE_DEFAULT_PARTITIONER == METIS_DEFAULT_PARTITIONER
969 ServiceProperty("DefaultPartitioner",ST_SubDomain),
972ARCANE_REGISTER_SERVICE_METISMESHPARTITIONER(DefaultPartitioner,MetisMeshPartitioner);
973#endif
974
975/*---------------------------------------------------------------------------*/
976/*---------------------------------------------------------------------------*/
977
978} // End namespace Arcane
979
980/*---------------------------------------------------------------------------*/
981/*---------------------------------------------------------------------------*/
#define ARCANE_THROW(exception_class,...)
Macro pour envoyer une exception avec formattage.
#define ARCANE_FATAL(...)
Macro envoyant une exception FatalErrorException.
#define ENUMERATE_(type, name, group)
Enumérateur générique d'un groupe d'entité
#define ENUMERATE_CELL(name, group)
Enumérateur générique d'un groupe de mailles.
#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.
ArcaneMetisMeshPartitionerObject(const Arcane::ServiceBuildInfo &sbi)
Constructeur.
CaseOptionsMetisMeshPartitioner * options() const
Options du jeu de données du service.
Exception lorsqu'un argument est invalide.
Conversion d'un tableau d'un type vers un autre type.
Vue modifiable d'un tableau d'un type T.
constexpr Integer size() const noexcept
Retourne la taille du tableau.
void fill(const DataType &data)
Remplissage du tableau.
ConstArrayView< T > constView() const
Vue constante sur ce 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.
const T * data() const
Accès à la racine du tableau hors toute protection.
iterator begin()
Itérateur sur le premier élément du tableau.
void add(ConstReferenceType val)
Ajoute l'élément val à la fin du tableau.
Maille d'un maillage.
Definition Item.h:1205
Int32 nbEdge() const
Nombre d'arêtes de la maille.
Definition Item.h:1295
Int32 nbFace() const
Nombre de faces de la maille.
Definition Item.h:1280
Int32 globalIteration() const
Numéro de l'itération courante.
Vue constante d'un tableau de type T.
Classe permettant d'activer/désactiver temporairement les exceptions flottantes du processeur.
Redistribute graph data to another "communicator".
void initWithOneRankPerNode(bool allow_only_one_rank)
Automatic distribution : one partitioning process per node.
Interface d'un partitionneur de maillage.
Interface d'un partitionneur de maillage.
Interface du gestionnaire de parallélisme pour un sous-domaine.
virtual bool isThreadImplementation() const =0
Indique si l'implémentation utilise les threads.
virtual void computeMinMaxSum(char val, char &min_val, char &max_val, char &sum_val, Int32 &min_rank, Int32 &max_rank)=0
Calcule en une opération la somme, le min, le max d'une valeur.
virtual Int32 commRank() const =0
Rang de cette instance dans le communicateur.
virtual Int32 commSize() const =0
Nombre d'instance dans le communicateur.
virtual void allGather(ConstArrayView< char > send_buf, ArrayView< char > recv_buf)=0
Effectue un regroupement sur tous les processeurs. Il s'agit d'une opération collective....
virtual Parallel::Communicator communicator() const =0
Communicateur MPI associé à ce gestionnaire.
virtual char reduce(eReduceType rt, char v)=0
Effectue la réduction de type rt sur le réel v et retourne la valeur.
virtual void barrier()=0
Effectue une barière.
Interface du gestionnaire d'un sous-domaine.
Definition ISubDomain.h:74
virtual const CommonVariables & commonVariables() const =0
Informations sur les variables standards.
virtual void flush()=0
Flush tous les flots.
@ PNoDump
Indique que la variable ne doit pas être sauvegardée.
Definition IVariable.h:74
bool hasFlags(Int32 flags) const
Retourne si les flags flags sont positionnées pour l'entité
Definition Item.h:338
Real imbalance() const override
Déséquilibre de temps de calcul.
Real maxImbalance() const override
Déséquilibre maximal autorisé
virtual void changeOwnersFromCells()
Positionne les nouveaux propriétaires des noeuds, arêtes et faces à partir des mailles.
Partitioneur de maillage utilisant la bibliothèque PARMetis.
void partitionMesh(bool initial_partition) override
void _removeEmptyPartsV2(Int32 nb_part, ArrayView< idxtype > metis_part)
Comble les partitions vides (version 2).
Int32 _removeEmptyPartsV2Helper(Int32 nb_part, ArrayView< idxtype > metis_part, Int32 algo_iteration)
Applique une itération de l'algorithme de suppression des partitions vides.
void build() override
Construction de niveau build du service.
void _partitionMesh(bool initial_partition, Int32 nb_part)
int _writeGraph(IParallelMng *pm, Span< const idxtype > metis_vtkdist, Span< const idxtype > metis_xadj, Span< const idxtype > metis_adjncy, Span< const idxtype > metis_vwgt, Span< const idxtype > metis_ewgt, const String &Name) const
void _removeEmptyPartsV1(Int32 nb_part, Int32 nb_own_cell, ArrayView< idxtype > metis_part)
Comble les partitions vides (version 1).
Wrapper autour des appels de Parmetis.
int callPartKway(const bool print_digest, const bool gather, idx_t *vtxdist, idx_t *xadj, idx_t *adjncy, idx_t *vwgt, idx_t *adjwgt, idx_t *wgtflag, idx_t *numflag, idx_t *ncon, idx_t *nparts, real_t *tpwgts, real_t *ubvec, idx_t *options, idx_t *edgecut, idx_t *part)
Simple wrapper autour de la routine ParMetis "ParMETIS_V3_PartKway".
int callAdaptiveRepart(const bool print_digest, const bool gather, idx_t *vtxdist, idx_t *xadj, idx_t *adjncy, idx_t *vwgt, idx_t *vsize, idx_t *adjwgt, idx_t *wgtflag, idx_t *numflag, idx_t *ncon, idx_t *nparts, real_t *tpwgts, real_t *ubvec, real_t *ipc2redist, idx_t *options, idx_t *edgecut, idx_t *part)
Simple wrapper autour de la routine ParMetis "ParMETIS_V3_AdaptiveRepart".
Conversion d'un tableau de flottants vers un tableau d'entiers/longs. \abstract Cette classe gere le ...
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.
constexpr __host__ __device__ SizeType size() const noexcept
Retourne la taille du tableau.
Definition Span.h:212
Vue d'un tableau d'éléments de type T.
Definition Span.h:513
Constructeur de chaîne de caractère unicode.
String toString() const
Retourne la chaîne de caractères construite.
Chaîne de caractères unicode.
bool null() const
Retourne true si la chaîne est nulle.
Definition String.cc:304
String upper() const
Transforme tous les caractères de la chaîne en majuscules.
Definition String.cc:466
const char * localstr() const
Retourne la conversion de l'instance dans l'encodage UTF-8.
Definition String.cc:227
Postionne le nom de l'action en cours d'exécution.
Definition Timer.h:110
TraceMessage info() const
Flot pour un message d'information.
ITraceMng * traceMng() const
Gestionnaire de trace.
TraceMessage pwarning() const
Vecteur 1D de données avec sémantique par valeur (style STL).
Paramètres nécessaires à la construction d'une variable.
ItemGroupT< Cell > CellGroup
Groupe de mailles.
Definition ItemTypes.h:183
#define ARCANE_REGISTER_SERVICE(aclass, a_service_property,...)
Macro pour enregistrer un service.
MeshVariableScalarRefT< Cell, Integer > VariableCellInteger
Grandeur au centre des mailles de type entier.
ItemVariableScalarRefT< Int32 > VariableItemInt32
Grandeur de type entier 32 bits.
Integer toInteger(Real r)
Converti un Int64 en un Integer.
Int32 toInt32(Int64 v)
Converti un Int64 en un Int32.
Integer toInteger(Real r)
Converti un Real en Integer.
Definition Convert.h:45
@ ReduceSum
Somme des valeurs.
@ ReduceMax
Maximum des valeurs.
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:426
std::int64_t Int64
Type entier signé sur 64 bits.
Int32 Integer
Type représentant un entier.
@ ST_SubDomain
Le service s'utilise au niveau du sous-domaine.
@ 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.