Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
MetisGraphGather.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/* MetisGraphGather.cc (C) 2000-2024 */
9/* */
10/* Graph gathering for 'Parmetis'. */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arcane/std/internal/MetisGraphGather.h"
15
16#include "arcane/utils/CheckedConvert.h"
17#include "arcane/utils/FatalErrorException.h"
18
19#include "arcane/core/IParallelMng.h"
20
21#include <algorithm>
22
23/*---------------------------------------------------------------------------*/
24/*---------------------------------------------------------------------------*/
25
26namespace
27{
28MPI_Comm _getMPICommunicator(Arcane::IParallelMng* pm)
29{
30 return static_cast<MPI_Comm>(pm->communicator());
31}
32} // namespace
33
34namespace Arcane
35{
36
37/*---------------------------------------------------------------------------*/
38/*---------------------------------------------------------------------------*/
39
40MetisGraphGather::
41MetisGraphGather(IParallelMng* pm)
42: TraceAccessor(pm->traceMng())
43, m_parallel_mng(pm)
44, m_my_rank(pm->commRank())
45, m_nb_rank(pm->commSize())
46{}
47
48/*---------------------------------------------------------------------------*/
49/*---------------------------------------------------------------------------*/
50
51/*---------------------------------------------------------------------------*/
52/*---------------------------------------------------------------------------*/
53
54template <class SourceType, class TargetType> void MetisGraphGather::
55_convertVector(const int size, ConstArrayView<SourceType> src, ArrayView<TargetType> dest)
56{
57 if (size > src.size())
58 ARCANE_FATAL("Source size is too small size={0} src_size={1}", size, src.size());
59 if (size > dest.size())
60 ARCANE_FATAL("Target size is too small size={0} target_size={1}", size, dest.size());
61 for (int i = 0; i < size; ++i) {
62 dest[i] = static_cast<TargetType>(src[i]);
63 }
64}
65
66/*---------------------------------------------------------------------------*/
67/*---------------------------------------------------------------------------*/
68
69void MetisGraphGather::
70gatherGraph(const bool need_part, ConstArrayView<idx_t> vtxdist, const int ncon, MetisGraphView my_graph,
71 MetisGraph& graph)
72{
73 MPI_Comm comm = _getMPICommunicator(m_parallel_mng);
74 info() << "Metis: gather graph";
75
76 const Int32 io_rank = m_parallel_mng->masterIORank();
77 const Int32 nb_rank = m_nb_rank;
78 const bool is_master_io = m_parallel_mng->isMasterIO();
79 // buffers
80
81 UniqueArray<int> my_buffer; // for the inputs of the MPI routines
82 UniqueArray<int> buffer; // for the outputs of the MPI routines
83 UniqueArray<int> offset; // for the gather
84
85 // number of vertices in the complete graph
86
87 if (m_my_rank == io_rank) {
88 graph.nb_vertices = CheckedConvert::toInt32(vtxdist[m_nb_rank]);
89 graph.have_vsize = my_graph.have_vsize; // it is assumed that all processors have the same value
90 graph.have_adjwgt = my_graph.have_adjwgt; // it is assumed that all processors have the same value
91 }
92 else {
93 graph.nb_vertices = 0;
94 graph.have_vsize = false;
95 graph.have_adjwgt = false;
96 }
97
98 // retrieves the dimensions characterizing the graph distribution across processors
99
100 my_buffer.resize(2);
101 if (m_my_rank == io_rank) {
102 offset.resize(nb_rank);
103 buffer.resize(2 * nb_rank);
104 }
105
106 my_buffer[0] = my_graph.nb_vertices; // number of vertices
107 my_buffer[1] = my_graph.adjncy.size(); // size of the adjacency list
108
109 MPI_Gather((void*)my_buffer.data(), 2, MPI_INT, (void*)buffer.data(), 2, MPI_INT, io_rank, comm);
110
111 int adjcny_size = 0;
112 UniqueArray<int> nb_vertices_per_rank(nb_rank);
113 UniqueArray<int> adjncy_size_per_rank(nb_rank);
114 UniqueArray<int> nb_vertice_cons_per_rank(nb_rank);
115 if (is_master_io) {
116 for (int rank = 0; rank < nb_rank; ++rank) {
117 nb_vertices_per_rank[rank] = buffer[2 * rank];
118 adjncy_size_per_rank[rank] = buffer[2 * rank + 1];
119 nb_vertice_cons_per_rank[rank] = nb_vertices_per_rank[rank] * ncon;
120 adjcny_size += adjncy_size_per_rank[rank];
121 }
122 }
123
124 // resize the buffers now that the dimensions are known (the largest array is adjcny)
125
126 int max_buffer_size = my_graph.nb_vertices * ncon;
127 max_buffer_size = std::max(my_graph.adjncy.size(), max_buffer_size);
128 my_buffer.resize(max_buffer_size);
129 if (is_master_io) {
130 max_buffer_size = graph.nb_vertices * ncon;
131 max_buffer_size = std::max(adjcny_size, max_buffer_size);
132 buffer.resize(max_buffer_size);
133 }
134
135 // Retrieves the adjacency list.
136
137 if (is_master_io) {
138 offset[0] = 0;
139 for (int rank = 1; rank < nb_rank; ++rank) {
140 offset[rank] = offset[rank - 1] + adjncy_size_per_rank[rank - 1];
141 }
142 graph.adjncy.resize(adjcny_size);
143 }
144
145 _convertVector(my_graph.adjncy.size(), my_graph.adjncy.constView(), my_buffer.view());
146
147 MPI_Gatherv(my_buffer.data(), my_graph.adjncy.size(), MPI_INT,
148 buffer.data(), adjncy_size_per_rank.data(),
149 offset.data(), MPI_INT, io_rank, comm);
150
151 _convertVector(adjcny_size, buffer.constView(), graph.adjncy.view());
152
153 // Retrieves the list of edge weights.
154
155 if (my_graph.have_adjwgt) {
156 if (is_master_io) {
157 graph.adjwgt.resize(adjcny_size);
158 }
159
160 _convertVector(my_graph.adjwgt.size(), my_graph.adjwgt.constView(), my_buffer.view());
161
162 MPI_Gatherv(my_buffer.data(), my_graph.adjwgt.size(), MPI_INT,
163 buffer.data(), adjncy_size_per_rank.data(),
164 offset.data(), MPI_INT, io_rank, comm);
165
166 _convertVector(adjcny_size, buffer.constView(), graph.adjwgt.view());
167 }
168
169 // Retrieves the list of indices in the adjacency list.
170 // At this stage, the retrieved list is not correct, except for the last index.
171 // In fact, the indices per processor are not the same as the indices
172 // in the complete adjacency list.
173
174 if (is_master_io) {
175 graph.xadj.resize(graph.nb_vertices + 1);
176 graph.xadj[graph.nb_vertices] = graph.adjncy.size();
177 offset[0] = 0;
178 for (int rank = 1; rank < nb_rank; ++rank) {
179 offset[rank] = offset[rank - 1] + nb_vertices_per_rank[rank - 1];
180 }
181 }
182
183 _convertVector(my_graph.nb_vertices, my_graph.xadj.constView(), my_buffer.view());
184
185 MPI_Gatherv(my_buffer.data(), my_graph.nb_vertices, MPI_INT,
186 buffer.data(), nb_vertices_per_rank.data(),
187 offset.data(), MPI_INT, io_rank, comm);
188
189 _convertVector(graph.nb_vertices, buffer.constView(), graph.xadj.view());
190
191 // Index correction
192
193 if (is_master_io) {
194 int start_adjncy_index = 0;
195 for (int rank = 1; rank < nb_rank; ++rank) {
196 start_adjncy_index += adjncy_size_per_rank[rank - 1];
197 //std::cerr << "rank " << rank << " offset " << start_adjncy_index << " vtxdist[rank] " << vtxdist[rank] << " vtxdist[rank+1] " << vtxdist[rank+1] << std::endl;
198 Int32 vtxdist_rank = CheckedConvert::toInt32(vtxdist[rank]);
199 Int32 vtxdist_rank_plus_one = CheckedConvert::toInt32(vtxdist[rank + 1]);
200 for (Int32 ixadj = vtxdist_rank; ixadj < vtxdist_rank_plus_one; ++ixadj) {
201 graph.xadj[ixadj] += start_adjncy_index;
202 }
203 }
204 }
205
206 // Retrieves the "memory" weights for the graph vertices.
207
208 if (my_graph.have_vsize) {
209 if (is_master_io) {
210 graph.vsize.resize(graph.nb_vertices);
211 }
212
213 _convertVector(my_graph.nb_vertices, my_graph.vsize.constView(), my_buffer.view());
214
215 MPI_Gatherv(my_buffer.data(), my_graph.nb_vertices, MPI_INT,
216 buffer.data(), nb_vertices_per_rank.data(),
217 offset.data(), MPI_INT, io_rank, comm);
218
219 _convertVector(graph.nb_vertices, buffer.constView(), graph.vsize.view());
220 }
221
222 // Retrieves the list of subdomain numbers from a previous partitioning
223
224 if (is_master_io) {
225 graph.part.resize(graph.nb_vertices);
226 }
227
228 if (need_part) {
229
230 _convertVector(my_graph.nb_vertices, my_graph.part.constView(), my_buffer.view());
231
232 MPI_Gatherv(my_buffer.data(), my_graph.nb_vertices, MPI_INT,
233 buffer.data(), nb_vertices_per_rank.data(),
234 offset.data(), MPI_INT, io_rank, comm);
235
236 _convertVector(graph.nb_vertices, buffer.constView(), graph.part.view());
237 }
238
239 // Retrieves the list of weights for the graph vertices. There may be several (ncon).
240
241 if (is_master_io) {
242 graph.vwgt.resize(graph.nb_vertices * ncon);
243 for (auto& x : offset) {
244 x *= ncon;
245 }
246 }
247
248 _convertVector(my_graph.nb_vertices * ncon, my_graph.vwgt.constView(), my_buffer.view());
249
250 MPI_Gatherv(my_buffer.data(), my_graph.nb_vertices * ncon, MPI_INT,
251 buffer.data(), nb_vertice_cons_per_rank.data(),
252 offset.data(), MPI_INT, io_rank, comm);
253
254 _convertVector(graph.nb_vertices * ncon, buffer.constView(), graph.vwgt.view());
255}
256
257/*---------------------------------------------------------------------------*/
258/*---------------------------------------------------------------------------*/
259
260void MetisGraphGather::
261scatterPart(ConstArrayView<idx_t> vtxdist, ConstArrayView<idx_t> part, ArrayView<idx_t> my_part)
262{
263 MPI_Comm comm = _getMPICommunicator(m_parallel_mng);
264
265 const Int32 nb_rank = m_nb_rank;
266 const bool is_master_io = m_parallel_mng->isMasterIO();
267 int io_rank = 0;
268
269 UniqueArray<Int32> nb_vertices(nb_rank);
270 UniqueArray<Int32> displ(nb_rank);
271
272 for (int rank = 0; rank < nb_rank; ++rank) {
273 displ[rank] = CheckedConvert::toInt32(vtxdist[rank]);
274 nb_vertices[rank] = CheckedConvert::toInt32(vtxdist[rank + 1] - vtxdist[rank]);
275 }
276
277 UniqueArray<int> send_buffer;
278 UniqueArray<int> recv_buffer(my_part.size());
279
280 if (is_master_io) {
281 send_buffer.resize(part.size());
282 }
283
284 _convertVector(send_buffer.size(), part, send_buffer.view());
285
286 MPI_Scatterv(send_buffer.data(), nb_vertices.data(),
287 displ.data(), MPI_INT, recv_buffer.data(),
288 my_part.size(), MPI_INT, io_rank, comm);
289
290 _convertVector(recv_buffer.size(), recv_buffer.constView(), my_part);
291}
292
293/*---------------------------------------------------------------------------*/
294/*---------------------------------------------------------------------------*/
295
296} // End namespace Arcane
297
298/*---------------------------------------------------------------------------*/
299/*---------------------------------------------------------------------------*/
#define ARCANE_FATAL(...)
Macro throwing a FatalErrorException.
Integer size() const
Number of elements in the vector.
Modifiable view of an array of type T.
constexpr Integer size() const noexcept
Returns the size of the array.
void resize(Int64 s)
Changes the number of elements in the array to s.
const T * data() const
Access to the root of the array without any protection.
ConstArrayView< T > constView() const
Constant view of this array.
ArrayView< T > view() const
Mutable view of this array.
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.
TraceMessage info() const
Flow for an information message.
1D data vector with value semantics (STL style).
Integer size() const
Number of elements in the vector.
constexpr ConstArrayView< T > constView() const noexcept
Constant view of this view.
constexpr Integer size() const noexcept
Returns the size of the array.
void resize(Int64 s)
Changes the number of elements in the array to s.
ArrayView< T > view() const
Mutable view of this array.
-- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature --
std::int32_t Int32
Signed integer type of 32 bits.