Arcane  v3.16.0.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
MetisWrapper.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2024 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/* MetisWrapper.cc (C) 2000-2024 */
9/* */
10/* Wrapper autour des appels de Parmetis. */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arcane/std/internal/MetisWrapper.h"
15
16#include "arcane/utils/CheckedConvert.h"
17#include "arcane/utils/ITraceMng.h"
18
19#include "arcane/core/IParallelMng.h"
20#include "arcane/core/ParallelMngUtils.h"
21
22#include "arcane/std/internal/MetisGraph.h"
23#include "arcane/std/internal/MetisGraphDigest.h"
24#include "arcane/std/internal/MetisGraphGather.h"
25
26#include <functional>
27
28/*---------------------------------------------------------------------------*/
29/*---------------------------------------------------------------------------*/
30
31namespace Arcane
32{
33
34/*---------------------------------------------------------------------------*/
35/*---------------------------------------------------------------------------*/
36
37MetisWrapper::
38MetisWrapper(IParallelMng* pm)
39: TraceAccessor(pm->traceMng())
40, m_parallel_mng(pm)
41{
42}
43
44/*---------------------------------------------------------------------------*/
45/*---------------------------------------------------------------------------*/
46
47/*---------------------------------------------------------------------------*/
48/*---------------------------------------------------------------------------*/
52int MetisWrapper::
53_callMetisWith2Processors(const Int32 ncon, const bool need_part,
54 ConstArrayView<idx_t> vtxdist, MetisGraphView my_graph,
55 MetisCall& metis_call)
56{
57 Int32 nb_rank = m_parallel_mng->commSize();
58 Int32 my_rank = m_parallel_mng->commRank();
59
60 UniqueArray<idx_t> half_vtxdist(vtxdist.size());
61
62 int comm_0_size = nb_rank / 2 + nb_rank % 2;
63 int comm_1_size = nb_rank / 2;
64 int comm_0_io_rank = 0;
65 int comm_1_io_rank = comm_0_size;
66
67 for (int i = 0; i < nb_rank + 1; ++i) {
68 half_vtxdist[i] = vtxdist[i];
69 }
70
71 int color = 1;
72 int key = my_rank;
73
74 if (my_rank >= comm_0_size) {
75 color = 2;
76 for (int i = 0; i < comm_1_size + 1; ++i) {
77 half_vtxdist[i] = vtxdist[i + comm_0_size] - vtxdist[comm_0_size];
78 }
79 }
80
81 MetisGraph metis_graph;
82 Ref<IParallelMng> half_pm = ParallelMngUtils::createSubParallelMngRef(m_parallel_mng, color, key);
83 MetisGraphGather metis_gather(half_pm.get());
84
85 metis_gather.gatherGraph(need_part, half_vtxdist, ncon, my_graph, metis_graph);
86
87 color = -1;
88 if (my_rank == comm_0_io_rank || my_rank == comm_1_io_rank) {
89 color = 1;
90 }
91
92 Ref<IParallelMng> metis_pm = ParallelMngUtils::createSubParallelMngRef(m_parallel_mng, color, key);
93
94 UniqueArray<idx_t> metis_vtxdist(3);
95 metis_vtxdist[0] = 0;
96 metis_vtxdist[1] = vtxdist[comm_0_size];
97 metis_vtxdist[2] = vtxdist[vtxdist.size() - 1];
98
99 int ierr = 0;
100
101 if (metis_pm.get()) {
102 MetisGraphView metis_graph_view(metis_graph);
103 ierr = metis_call(metis_pm.get(), metis_graph_view, metis_vtxdist);
104 }
105
106 // S'assure que tout le monde a la même valeur de l'erreur.
107 half_pm->broadcast(ArrayView<int>(1, &ierr), 0);
108
109 metis_gather.scatterPart(half_vtxdist, metis_graph.part, my_graph.part);
110
111 return ierr;
112}
113
114/*---------------------------------------------------------------------------*/
115/*---------------------------------------------------------------------------*/
116
117int MetisWrapper::
118callPartKway(const bool print_digest, const bool gather,
119 idx_t *vtxdist, idx_t *xadj, idx_t *adjncy, idx_t *vwgt,
120 idx_t *adjwgt, idx_t *wgtflag, idx_t *numflag, idx_t *ncon, idx_t *nparts,
121 real_t *tpwgts, real_t *ubvec, idx_t *options, idx_t *edgecut, idx_t *part)
122{
123 int ierr = 0;
124 Int32 nb_rank = m_parallel_mng->commSize();
125 Int32 my_rank = m_parallel_mng->commRank();
126
127 MetisCall partkway = [&](IParallelMng* pm, MetisGraphView graph,
128 ArrayView<idx_t> graph_vtxdist)
129 {
130 MPI_Comm graph_comm = static_cast<MPI_Comm>(pm->communicator());
131 // NOTE GG: il peut arriver que ces deux pointeurs soient nuls s'il n'y a pas
132 // de voisins. Si tout le reste est cohérent cela ne pose pas de problèmes mais ParMetis
133 // n'aime pas les tableaux vides donc si c'est le cas on met un 0.
134 // TODO: il faudrait regarder en amont s'il ne vaudrait pas mieux mettre des valeurs
135 // dans ces deux tableaux au cas où.
136 idx_t null_idx = 0;
137 idx_t* adjncy_data = graph.adjncy.data();
138 idx_t* adjwgt_data = graph.adjwgt.data();
139 if (!adjncy_data)
140 adjncy_data = &null_idx;
141 if (!adjwgt_data)
142 adjwgt_data = &null_idx;
143 return ParMETIS_V3_PartKway(graph_vtxdist.data(), graph.xadj.data(),
144 adjncy_data, graph.vwgt.data(),
145 adjwgt_data, wgtflag, numflag, ncon, nparts, tpwgts,
146 ubvec, options, edgecut, graph.part.data(), &graph_comm);
147 };
148
149 // Version séquentielle utilisant directement Metis
150 // NOTE: Ce bout de code est le même que dans callAdaptativeRepart.
151 // Il faudrait mutualiser les deux.
152 MetisCall partkway_seq = [&](IParallelMng*, MetisGraphView graph,
153 ArrayView<idx_t> graph_vtxdist)
154 {
155 idx_t options2[METIS_NOPTIONS];
156 METIS_SetDefaultOptions(options2);
157 options2[METIS_OPTION_CTYPE] = METIS_CTYPE_SHEM;
158 options2[METIS_OPTION_UFACTOR] = 30; // TODO: a récupérer depuis les options de parmetis
159 options2[METIS_OPTION_NUMBERING] = 0; // 0 pour indiquer que les tableaux sont en C (indices commencent à 0)
160 options2[METIS_OPTION_MINCONN] = 0;
161 options2[METIS_OPTION_OBJTYPE] = METIS_OBJTYPE_VOL;
162 options2[METIS_OPTION_SEED] = 25; // TODO: pouvoir changer la graine
163 pwarning() << "MetisWrapper: using user 'imbalance_factor' is not yet implemented. Using defaut value 30";
164
165 // Le nombre de sommets du graph est dans le premier indice de graph_vtxdist
166 idx_t nvtxs = graph_vtxdist[1];
167 return METIS_PartGraphKway(&nvtxs /*graph_vtxdist.data()*/, ncon, graph.xadj.data(),
168 graph.adjncy.data(), graph.vwgt.data(), graph.vsize.data(),
169 graph.adjwgt.data(), nparts, tpwgts,
170 ubvec, options2, edgecut, graph.part.data());
171 };
172
173 MetisGraphView my_graph;
174
175 ArrayView<idx_t> offset(nb_rank + 1, vtxdist);
176 my_graph.nb_vertices = CheckedConvert::toInt32(offset[my_rank + 1] - offset[my_rank]);
177 my_graph.xadj = ArrayView<idx_t>(my_graph.nb_vertices + 1, xadj);
178 const Int32 adjacency_size = CheckedConvert::toInt32(my_graph.xadj[my_graph.nb_vertices]);
179 const Int32 nb_con = CheckedConvert::toInt32(*ncon);
180 my_graph.adjncy = ArrayView<idx_t>(adjacency_size, adjncy);
181 my_graph.vwgt = ArrayView<idx_t>(CheckedConvert::multiply(my_graph.nb_vertices, nb_con), vwgt);
182 my_graph.adjwgt = ArrayView<idx_t>(adjacency_size, adjwgt);
183 my_graph.part = ArrayView<idx_t>(my_graph.nb_vertices, part);
184 my_graph.have_vsize = false;
185 my_graph.have_adjwgt = true;
186
187 if (print_digest){
188 MetisGraphDigest d(m_parallel_mng);
189 String digest = d.computeInputDigest(false, 3, my_graph, vtxdist, wgtflag, numflag,
190 ncon, nparts, tpwgts, ubvec, nullptr, options);
191 if (my_rank == 0) {
192 info() << "signature des entrees Metis = " << digest;
193 }
194 }
195
196 if (gather && nb_rank > 2) {
197 // Normalement c'est plus rapide ...
198 info() << "Partitioning metis : re-grouping " << nb_rank << " -> 2 rank";
199 ierr = _callMetisWith2Processors(nb_con, false, offset, my_graph, partkway);
200 }
201 else {
202 info() << "Partitioning metis : nb rank = " << nb_rank;
203 MetisCall& metis_call = (nb_rank == 1) ? partkway_seq : partkway;
204 ierr = metis_call(m_parallel_mng, my_graph, offset);
205 }
206
207 info() << "End Partitioning metis";
208 if (print_digest){
209 MetisGraphDigest d(m_parallel_mng);
210 String digest = d.computeOutputDigest(my_graph, edgecut);
211 if (my_rank == 0) {
212 info() << "hash for Metis output = " << digest;
213 }
214 }
215
216 return ierr;
217}
218
219/*---------------------------------------------------------------------------*/
220/*---------------------------------------------------------------------------*/
221
222int MetisWrapper::
223callAdaptiveRepart(const bool print_digest, const bool gather,
224 idx_t *vtxdist, idx_t *xadj, idx_t *adjncy, idx_t *vwgt,
225 idx_t *vsize, idx_t *adjwgt, idx_t *wgtflag, idx_t *numflag, idx_t *ncon,
226 idx_t *nparts, real_t *tpwgts, real_t *ubvec, real_t *ipc2redist,
227 idx_t *options, idx_t *edgecut, idx_t *part)
228{
229 int ierr = 0;
230 Int32 nb_rank = m_parallel_mng->commSize();
231 Int32 my_rank = m_parallel_mng->commRank();
232
233 MetisCall repart_func = [&](IParallelMng* pm, MetisGraphView graph,
234 ArrayView<idx_t> graph_vtxdist)
235 {
236 MPI_Comm graph_comm = static_cast<MPI_Comm>(pm->communicator());
237 return ParMETIS_V3_AdaptiveRepart(graph_vtxdist.data(), graph.xadj.data(),
238 graph.adjncy.data(), graph.vwgt.data(),
239 graph.vsize.data(), graph.adjwgt.data(),
240 wgtflag, numflag, ncon, nparts, tpwgts, ubvec,
241 ipc2redist, options, edgecut,
242 graph.part.data(), &graph_comm);
243 };
244
245 // Version séquentielle utilisant directement Metis
246 // NOTE: Ce bout de code est le même que dans callPartKWay
247 // Il faudrait mutualiser les deux.
248 MetisCall repart_seq_func = [&](IParallelMng*, MetisGraphView graph,
249 ArrayView<idx_t> graph_vtxdist)
250 {
251 idx_t options2[METIS_NOPTIONS];
252 METIS_SetDefaultOptions(options2);
253 options2[METIS_OPTION_CTYPE] = METIS_CTYPE_SHEM;
254 options2[METIS_OPTION_UFACTOR] = 30; // TODO: a récupérer depuis les options de parmetis
255 options2[METIS_OPTION_NUMBERING] = 0; // 0 pour indiquer que les tableaux sont en C (indices commencent à 0)
256 options2[METIS_OPTION_MINCONN] = 0;
257 options2[METIS_OPTION_OBJTYPE] = METIS_OBJTYPE_VOL;
258 options2[METIS_OPTION_SEED] = 25; // TODO: pouvoir changer la graine
259 pwarning() << "MetisWrapper: using user 'imbalance_factor' is not yet implemented. Using defaut value 30";
260 // Le nombre de sommets du graph est dans le premier indice de graph_vtxdist
261 idx_t nvtxs = graph_vtxdist[1];
262 return METIS_PartGraphKway(&nvtxs /*graph_vtxdist.data()*/, ncon, graph.xadj.data(),
263 graph.adjncy.data(), graph.vwgt.data(), graph.vsize.data(),
264 graph.adjwgt.data(), nparts, tpwgts,
265 ubvec, options2, edgecut, graph.part.data());
266 };
267
268 MetisGraphView my_graph;
269
270
271 ArrayView<idx_t> offset(nb_rank + 1, vtxdist);
272 my_graph.nb_vertices = CheckedConvert::toInt32(offset[my_rank + 1] - offset[my_rank]);
273 my_graph.xadj = ArrayView<idx_t>(my_graph.nb_vertices + 1, xadj);
274 const Int32 adjacency_size = CheckedConvert::toInt32(my_graph.xadj[my_graph.nb_vertices]);
275 const Int32 nb_con = CheckedConvert::toInt32(*ncon);
276 my_graph.adjncy = ArrayView<idx_t>(adjacency_size, adjncy);
277 my_graph.vwgt = ArrayView<idx_t>(CheckedConvert::multiply(my_graph.nb_vertices, nb_con), vwgt);
278 my_graph.vsize = ArrayView<idx_t>(my_graph.nb_vertices, vsize);
279 my_graph.adjwgt = ArrayView<idx_t>(adjacency_size, adjwgt);
280 my_graph.part = ArrayView<idx_t>(my_graph.nb_vertices, part);
281 my_graph.have_vsize = true;
282 my_graph.have_adjwgt = true;
283
284
285 if (print_digest){
286 MetisGraphDigest d(m_parallel_mng);
287 String digest = d.computeInputDigest(true, 4, my_graph, vtxdist, wgtflag, numflag,
288 ncon, nparts, tpwgts, ubvec, nullptr, options);
289 if (my_rank == 0) {
290 info() << "signature des entrees Metis = " << digest;
291 }
292 }
293
294 if (gather && nb_rank > 2) {
295 info() << "Partionnement metis : regroupement " << nb_rank << " -> 2 processeurs";
296 ierr = _callMetisWith2Processors(nb_con, true, offset, my_graph, repart_func);
297 }
298 else {
299 info() << "Partionnement metis : nb processeurs = " << nb_rank;
300 MetisCall& metis_call = (nb_rank == 1) ? repart_seq_func : repart_func;
301 ierr = metis_call(m_parallel_mng, my_graph, offset);
302 }
303
304 if (print_digest) {
305 MetisGraphDigest d(m_parallel_mng);
306 String digest = d.computeOutputDigest(my_graph, edgecut);
307 if (my_rank == 0) {
308 info() << "signature des sorties Metis = " << digest;
309 }
310 }
311
312 return ierr;
313}
314
315/*---------------------------------------------------------------------------*/
316/*---------------------------------------------------------------------------*/
317
318} // End namespace Arcane
319
320/*---------------------------------------------------------------------------*/
321/*---------------------------------------------------------------------------*/
Vue modifiable d'un tableau d'un type T.
Vue constante d'un tableau de type T.
constexpr Integer size() const noexcept
Nombre d'éléments du tableau.
Interface du gestionnaire de parallélisme pour un sous-domaine.
virtual Parallel::Communicator communicator() const =0
Communicateur MPI associé à ce gestionnaire.
Calcule une somme de contrôle globale des entrées/sorties Metis.
String computeInputDigest(const bool need_part, const int nb_options, const MetisGraphView &my_graph, const idx_t *vtxdist, const idx_t *wgtflag, const idx_t *numflag, const idx_t *ncon, const idx_t *nparts, const real_t *tpwgts, const real_t *ubvec, const real_t *ipc2redist, const idx_t *options)
String computeOutputDigest(const MetisGraphView &my_graph, const idx_t *edgecut)
void scatterPart(ConstArrayView< idx_t > vtxdist, ConstArrayView< idx_t > part, ArrayView< idx_t > my_part)
Distribue le partitionnement "part" depuis le processeur de rang 0 dans le communicateur "comm" sur t...
void gatherGraph(const bool need_part, ConstArrayView< idx_t > vtxdist, const int ncon, MetisGraphView my_graph, MetisGraph &graph)
Effectue un regroupement du graphe ParMetis "my_graph" sur le processeur de rang 0 dans le communicat...
int _callMetisWith2Processors(const Int32 ncon, const bool need_part, ConstArrayView< idx_t > vtxdist, MetisGraphView my_graph, MetisCall &metis)
Appelle Metis en regroupant le graph sur 2 processeurs.
Référence à une instance.
InstanceType * get() const
Instance associée ou nullptr si aucune.
Chaîne de caractères unicode.
TraceMessage info() const
Flot pour un message d'information.
TraceMessage pwarning() const
Vecteur 1D de données avec sémantique par valeur (style STL).
Integer multiply(Integer x, Integer y, Integer z)
Multiplie trois 'Integer' et vérifie que le résultat peut être contenu dans un 'Integer'.
Int32 toInt32(Int64 v)
Converti un Int64 en un Int32.
Ref< IParallelMng > createSubParallelMngRef(IParallelMng *pm, Int32 color, Int32 key)
Créé un nouveau gestionnaire de parallélisme pour un sous-ensemble des rangs.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
std::int32_t Int32
Type entier signé sur 32 bits.