Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
MetisWrapper.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* MetisWrapper.cc (C) 2000-2024 */
9/* */
10/* Wrapper around Parmetis calls. */
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/*---------------------------------------------------------------------------*/
49
53int MetisWrapper::
54_callMetisWith2Processors(const Int32 ncon, const bool need_part,
55 ConstArrayView<idx_t> vtxdist, MetisGraphView my_graph,
56 MetisCall& metis_call)
57{
58 Int32 nb_rank = m_parallel_mng->commSize();
59 Int32 my_rank = m_parallel_mng->commRank();
60
61 UniqueArray<idx_t> half_vtxdist(vtxdist.size());
62
63 int comm_0_size = nb_rank / 2 + nb_rank % 2;
64 int comm_1_size = nb_rank / 2;
65 int comm_0_io_rank = 0;
66 int comm_1_io_rank = comm_0_size;
67
68 for (int i = 0; i < nb_rank + 1; ++i) {
69 half_vtxdist[i] = vtxdist[i];
70 }
71
72 int color = 1;
73 int key = my_rank;
74
75 if (my_rank >= comm_0_size) {
76 color = 2;
77 for (int i = 0; i < comm_1_size + 1; ++i) {
78 half_vtxdist[i] = vtxdist[i + comm_0_size] - vtxdist[comm_0_size];
79 }
80 }
81
82 MetisGraph metis_graph;
83 Ref<IParallelMng> half_pm = ParallelMngUtils::createSubParallelMngRef(m_parallel_mng, color, key);
84 MetisGraphGather metis_gather(half_pm.get());
85
86 metis_gather.gatherGraph(need_part, half_vtxdist, ncon, my_graph, metis_graph);
87
88 color = -1;
89 if (my_rank == comm_0_io_rank || my_rank == comm_1_io_rank) {
90 color = 1;
91 }
92
93 Ref<IParallelMng> metis_pm = ParallelMngUtils::createSubParallelMngRef(m_parallel_mng, color, key);
94
95 UniqueArray<idx_t> metis_vtxdist(3);
96 metis_vtxdist[0] = 0;
97 metis_vtxdist[1] = vtxdist[comm_0_size];
98 metis_vtxdist[2] = vtxdist[vtxdist.size() - 1];
99
100 int ierr = 0;
101
102 if (metis_pm.get()) {
103 MetisGraphView metis_graph_view(metis_graph);
104 ierr = metis_call(metis_pm.get(), metis_graph_view, metis_vtxdist);
105 }
106
107 // Ensures everyone has the same error value.
108 half_pm->broadcast(ArrayView<int>(1, &ierr), 0);
109
110 metis_gather.scatterPart(half_vtxdist, metis_graph.part, my_graph.part);
111
112 return ierr;
113}
114
115/*---------------------------------------------------------------------------*/
116/*---------------------------------------------------------------------------*/
117
118int MetisWrapper::
119callPartKway(const bool print_digest, const bool gather,
120 idx_t* vtxdist, idx_t* xadj, idx_t* adjncy, idx_t* vwgt,
121 idx_t* adjwgt, idx_t* wgtflag, idx_t* numflag, idx_t* ncon, idx_t* nparts,
122 real_t* tpwgts, real_t* ubvec, idx_t* options, idx_t* edgecut, idx_t* part)
123{
124 int ierr = 0;
125 Int32 nb_rank = m_parallel_mng->commSize();
126 Int32 my_rank = m_parallel_mng->commRank();
127
128 MetisCall partkway = [&](IParallelMng* pm, MetisGraphView graph,
129 ArrayView<idx_t> graph_vtxdist) {
130 MPI_Comm graph_comm = static_cast<MPI_Comm>(pm->communicator());
131 // NOTE GG: it may happen that these two pointers are null if there are no
132 // neighbors. If everything else is consistent, this is not a problem, but ParMetis
133 // does not like empty arrays, so if that is the case, we put a 0.
134 // TODO: we should check upstream if it wouldn't be better to put values
135 // in these two arrays in case.
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 // Sequential version using Metis directly
150 // NOTE: This piece of code is the same as in callAdaptiveRepart.
151 // We should centralize/abstract the two.
152 MetisCall partkway_seq = [&](IParallelMng*, MetisGraphView graph,
153 ArrayView<idx_t> graph_vtxdist) {
154 idx_t options2[METIS_NOPTIONS];
155 METIS_SetDefaultOptions(options2);
156 options2[METIS_OPTION_CTYPE] = METIS_CTYPE_SHEM;
157 options2[METIS_OPTION_UFACTOR] = 30; // TODO: retrieve from parmetis options
158 options2[METIS_OPTION_NUMBERING] = 0; // 0 to indicate that the arrays are in C (indices start at 0)
159 options2[METIS_OPTION_MINCONN] = 0;
160 options2[METIS_OPTION_OBJTYPE] = METIS_OBJTYPE_VOL;
161 options2[METIS_OPTION_SEED] = 25; // TODO: be able to change the seed
162 pwarning() << "MetisWrapper: using user 'imbalance_factor' is not yet implemented. Using defaut value 30";
163
164 // The number of vertices of the graph is in the first index of graph_vtxdist
165 idx_t nvtxs = graph_vtxdist[1];
166 return METIS_PartGraphKway(&nvtxs /*graph_vtxdist.data()*/, ncon, graph.xadj.data(),
167 graph.adjncy.data(), graph.vwgt.data(), graph.vsize.data(),
168 graph.adjwgt.data(), nparts, tpwgts,
169 ubvec, options2, edgecut, graph.part.data());
170 };
171
172 MetisGraphView my_graph;
173
174 ArrayView<idx_t> offset(nb_rank + 1, vtxdist);
175 my_graph.nb_vertices = CheckedConvert::toInt32(offset[my_rank + 1] - offset[my_rank]);
176 my_graph.xadj = ArrayView<idx_t>(my_graph.nb_vertices + 1, xadj);
177 const Int32 adjacency_size = CheckedConvert::toInt32(my_graph.xadj[my_graph.nb_vertices]);
178 const Int32 nb_con = CheckedConvert::toInt32(*ncon);
179 my_graph.adjncy = ArrayView<idx_t>(adjacency_size, adjncy);
180 my_graph.vwgt = ArrayView<idx_t>(CheckedConvert::multiply(my_graph.nb_vertices, nb_con), vwgt);
181 my_graph.adjwgt = ArrayView<idx_t>(adjacency_size, adjwgt);
182 my_graph.part = ArrayView<idx_t>(my_graph.nb_vertices, part);
183 my_graph.have_vsize = false;
184 my_graph.have_adjwgt = true;
185
186 if (print_digest) {
187 MetisGraphDigest d(m_parallel_mng);
188 String digest = d.computeInputDigest(false, 3, my_graph, vtxdist, wgtflag, numflag,
189 ncon, nparts, tpwgts, ubvec, nullptr, options);
190 if (my_rank == 0) {
191 info() << "Metis input signature = " << digest;
192 }
193 }
194
195 if (gather && nb_rank > 2) {
196 // Normally this is faster ...
197 info() << "Partitioning metis: re-grouping " << nb_rank << " -> 2 ranks";
198 ierr = _callMetisWith2Processors(nb_con, false, offset, my_graph, partkway);
199 }
200 else {
201 info() << "Partitioning metis: nb ranks = " << nb_rank;
202 MetisCall& metis_call = (nb_rank == 1) ? partkway_seq : partkway;
203 ierr = metis_call(m_parallel_mng, my_graph, offset);
204 }
205
206 info() << "End Partitioning metis";
207 if (print_digest) {
208 MetisGraphDigest d(m_parallel_mng);
209 String digest = d.computeOutputDigest(my_graph, edgecut);
210 if (my_rank == 0) {
211 info() << "hash for Metis output = " << digest;
212 }
213 }
214
215 return ierr;
216}
217
218/*---------------------------------------------------------------------------*/
219/*---------------------------------------------------------------------------*/
220
221int MetisWrapper::
222callAdaptiveRepart(const bool print_digest, const bool gather,
223 idx_t* vtxdist, idx_t* xadj, idx_t* adjncy, idx_t* vwgt,
224 idx_t* vsize, idx_t* adjwgt, idx_t* wgtflag, idx_t* numflag, idx_t* ncon,
225 idx_t* nparts, real_t* tpwgts, real_t* ubvec, real_t* ipc2redist,
226 idx_t* options, idx_t* edgecut, idx_t* part)
227{
228 int ierr = 0;
229 Int32 nb_rank = m_parallel_mng->commSize();
230 Int32 my_rank = m_parallel_mng->commRank();
231
232 MetisCall repart_func = [&](IParallelMng* pm, MetisGraphView graph,
233 ArrayView<idx_t> graph_vtxdist) {
234 MPI_Comm graph_comm = static_cast<MPI_Comm>(pm->communicator());
235 return ParMETIS_V3_AdaptiveRepart(graph_vtxdist.data(), graph.xadj.data(),
236 graph.adjncy.data(), graph.vwgt.data(),
237 graph.vsize.data(), graph.adjwgt.data(),
238 wgtflag, numflag, ncon, nparts, tpwgts, ubvec,
239 ipc2redist, options, edgecut,
240 graph.part.data(), &graph_comm);
241 };
242
243 // Sequential version using Metis directly
244 // NOTE: This piece of code is the same as in callPartKWay
245 // We should centralize/abstract the two.
246 MetisCall repart_seq_func = [&](IParallelMng*, MetisGraphView graph,
247 ArrayView<idx_t> graph_vtxdist) {
248 idx_t options2[METIS_NOPTIONS];
249 METIS_SetDefaultOptions(options2);
250 options2[METIS_OPTION_CTYPE] = METIS_CTYPE_SHEM;
251 options2[METIS_OPTION_UFACTOR] = 30; // TODO: retrieve from parmetis options
252 options2[METIS_OPTION_NUMBERING] = 0; // 0 to indicate that the arrays are in C (indices start at 0)
253 options2[METIS_OPTION_MINCONN] = 0;
254 options2[METIS_OPTION_OBJTYPE] = METIS_OBJTYPE_VOL;
255 options2[METIS_OPTION_SEED] = 25; // TODO: be able to change the seed
256 pwarning() << "MetisWrapper: using user 'imbalance_factor' is not yet implemented. Using defaut value 30";
257 // The number of vertices of the graph is in the first index of graph_vtxdist
258 idx_t nvtxs = graph_vtxdist[1];
259 return METIS_PartGraphKway(&nvtxs /*graph_vtxdist.data()*/, ncon, graph.xadj.data(),
260 graph.adjncy.data(), graph.vwgt.data(), graph.vsize.data(),
261 graph.adjwgt.data(), nparts, tpwgts,
262 ubvec, options2, edgecut, graph.part.data());
263 };
264
265 MetisGraphView my_graph;
266
267 ArrayView<idx_t> offset(nb_rank + 1, vtxdist);
268 my_graph.nb_vertices = CheckedConvert::toInt32(offset[my_rank + 1] - offset[my_rank]);
269 my_graph.xadj = ArrayView<idx_t>(my_graph.nb_vertices + 1, xadj);
270 const Int32 adjacency_size = CheckedConvert::toInt32(my_graph.xadj[my_graph.nb_vertices]);
271 const Int32 nb_con = CheckedConvert::toInt32(*ncon);
272 my_graph.adjncy = ArrayView<idx_t>(adjacency_size, adjncy);
273 my_graph.vwgt = ArrayView<idx_t>(CheckedConvert::multiply(my_graph.nb_vertices, nb_con), vwgt);
274 my_graph.vsize = ArrayView<idx_t>(my_graph.nb_vertices, vsize);
275 my_graph.adjwgt = ArrayView<idx_t>(adjacency_size, adjwgt);
276 my_graph.part = ArrayView<idx_t>(my_graph.nb_vertices, part);
277 my_graph.have_vsize = true;
278 my_graph.have_adjwgt = true;
279
280 if (print_digest) {
281 MetisGraphDigest d(m_parallel_mng);
282 String digest = d.computeInputDigest(true, 4, my_graph, vtxdist, wgtflag, numflag,
283 ncon, nparts, tpwgts, ubvec, nullptr, options);
284 if (my_rank == 0) {
285 info() << "Metis input signature = " << digest;
286 }
287 }
288
289 if (gather && nb_rank > 2) {
290 info() << "Metis partitioning: regrouping " << nb_rank << " -> 2 processors";
291 ierr = _callMetisWith2Processors(nb_con, true, offset, my_graph, repart_func);
292 }
293 else {
294 info() << "Metis partitioning: nb processors = " << nb_rank;
295 MetisCall& metis_call = (nb_rank == 1) ? repart_seq_func : repart_func;
296 ierr = metis_call(m_parallel_mng, my_graph, offset);
297 }
298
299 if (print_digest) {
300 MetisGraphDigest d(m_parallel_mng);
301 String digest = d.computeOutputDigest(my_graph, edgecut);
302 if (my_rank == 0) {
303 info() << "Metis output signature = " << digest;
304 }
305 }
306
307 return ierr;
308}
309
310/*---------------------------------------------------------------------------*/
311/*---------------------------------------------------------------------------*/
312
313} // End namespace Arcane
314
315/*---------------------------------------------------------------------------*/
316/*---------------------------------------------------------------------------*/
Modifiable view of an array of type T.
Constant view of an array of type T.
constexpr Integer size() const noexcept
Number of elements in the array.
Interface of the parallelism manager for a subdomain.
virtual Parallel::Communicator communicator() const =0
MPI communicator associated with this manager.
Calculates a global checksum of Metis inputs/outputs.
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)
Distributes the partitioning "part" from processor rank 0 in the communicator "comm" to all processor...
void gatherGraph(const bool need_part, ConstArrayView< idx_t > vtxdist, const int ncon, MetisGraphView my_graph, MetisGraph &graph)
Performs a gathering of the ParMetis graph "my_graph" on processor rank 0 in the communicator "comm"....
int _callMetisWith2Processors(const Int32 ncon, const bool need_part, ConstArrayView< idx_t > vtxdist, MetisGraphView my_graph, MetisCall &metis)
Calls Metis by grouping the graph onto 2 processors.
InstanceType * get() const
Associated instance or nullptr if none.
Reference to an instance.
TraceMessage info() const
Flow for an information message.
TraceMessage pwarning() const
1D data vector with value semantics (STL style).
Ref< IParallelMng > createSubParallelMngRef(IParallelMng *pm, Int32 color, Int32 key)
Creates a new parallelism manager for a subset of ranks.
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --
std::int32_t Int32
Signed integer type of 32 bits.