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