Arcane  v3.14.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
MetisGraphGather.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/* MetisGraphGather.cc (C) 2000-2024 */
9/* */
10/* Regroupement de graphes de '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,
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; // pour les entrees des routines MPI
82 UniqueArray<int> buffer; // pour les sorties des routines MPI
83 UniqueArray<int> offset; // pour les gather
84
85 // nombre de sommets du graph complet
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; // on suppose que tous les processeurs ont la meme valeur
90 graph.have_adjwgt = my_graph.have_adjwgt; // on suppose que tous les processeurs ont la meme valeur
91 } else {
92 graph.nb_vertices = 0;
93 graph.have_vsize = false;
94 graph.have_adjwgt = false;
95 }
96
97 // récupère les dimensions caractérisant la répartition du graphe sur les processeurs
98
99 my_buffer.resize(2);
100 if (m_my_rank == io_rank) {
101 offset.resize(nb_rank);
102 buffer.resize(2 * nb_rank);
103 }
104
105 my_buffer[0] = my_graph.nb_vertices; // nb de sommets
106 my_buffer[1] = my_graph.adjncy.size(); // taille de la liste d'adjacences
107
108 MPI_Gather((void*)my_buffer.data(), 2, MPI_INT, (void*)buffer.data(), 2, MPI_INT, io_rank, comm);
109
110 int adjcny_size = 0;
114 if (is_master_io) {
115 for (int rank = 0; rank < nb_rank; ++rank) {
116 nb_vertices_per_rank[rank] = buffer[2 * rank];
117 adjncy_size_per_rank[rank] = buffer[2 * rank + 1];
120 }
121 }
122
123 // retaille les buffers maintenant que les dimensions sont connues (le plus gros tableau est adjcny)
124
125 int max_buffer_size = my_graph.nb_vertices * ncon;
126 max_buffer_size = std::max(my_graph.adjncy.size(), max_buffer_size);
128 if (is_master_io) {
129 max_buffer_size = graph.nb_vertices * ncon;
131 buffer.resize(max_buffer_size);
132 }
133
134 // Récupère la liste d'adjacences.
135
136 if (is_master_io) {
137 offset[0] = 0;
138 for (int rank = 1; rank < nb_rank; ++rank) {
139 offset[rank] = offset[rank - 1] + adjncy_size_per_rank[rank - 1];
140 }
141 graph.adjncy.resize(adjcny_size);
142 }
143
144 _convertVector(my_graph.adjncy.size(), my_graph.adjncy.constView(), my_buffer.view());
145
146 MPI_Gatherv(my_buffer.data(), my_graph.adjncy.size(), MPI_INT,
147 buffer.data(), adjncy_size_per_rank.data(),
148 offset.data(), MPI_INT, io_rank, comm);
149
150 _convertVector(adjcny_size, buffer.constView(), graph.adjncy.view());
151
152 // Recupere la liste des poids aux arretes du graph.
153
154 if (my_graph.have_adjwgt) {
155 if (is_master_io) {
156 graph.adjwgt.resize(adjcny_size);
157 }
158
159 _convertVector(my_graph.adjwgt.size(), my_graph.adjwgt.constView(), my_buffer.view());
160
161 MPI_Gatherv(my_buffer.data(), my_graph.adjwgt.size(), MPI_INT,
162 buffer.data(), adjncy_size_per_rank.data(),
163 offset.data(), MPI_INT, io_rank, comm);
164
165 _convertVector(adjcny_size, buffer.constView(), graph.adjwgt.view());
166 }
167
168 // Récupère la liste des index dans la liste d'adjacences.
169 // A cette étape, la liste récupérée n'est pas correcte, sauf le dernier indice.
170 // En effet les indexes par processeur ne sont pas les meme que les index
171 // dans la liste d'adjacences complete.
172
173 if (is_master_io) {
174 graph.xadj.resize(graph.nb_vertices + 1);
175 graph.xadj[graph.nb_vertices] = graph.adjncy.size();
176 offset[0] = 0;
177 for (int rank = 1; rank < nb_rank; ++rank) {
178 offset[rank] = offset[rank - 1] + nb_vertices_per_rank[rank - 1];
179 }
180 }
181
182 _convertVector(my_graph.nb_vertices, my_graph.xadj.constView(), my_buffer.view());
183
184 MPI_Gatherv(my_buffer.data(), my_graph.nb_vertices, MPI_INT,
185 buffer.data(), nb_vertices_per_rank.data(),
186 offset.data(), MPI_INT, io_rank, comm);
187
188 _convertVector(graph.nb_vertices, buffer.constView(), graph.xadj.view());
189
190 // Correction des index
191
192 if (is_master_io) {
193 int start_adjncy_index = 0;
194 for (int rank = 1; rank < nb_rank; ++rank) {
196 //std::cerr << "rank " << rank << " offset " << start_adjncy_index << " vtxdist[rank] " << vtxdist[rank] << " vtxdist[rank+1] " << vtxdist[rank+1] << std::endl;
197 Int32 vtxdist_rank = CheckedConvert::toInt32(vtxdist[rank]);
198 Int32 vtxdist_rank_plus_one = CheckedConvert::toInt32(vtxdist[rank + 1]);
201 }
202 }
203 }
204
205 // Récupère la liste des poids "memoire" aux sommets du graph.
206
207 if (my_graph.have_vsize) {
208 if (is_master_io) {
209 graph.vsize.resize(graph.nb_vertices);
210 }
211
212 _convertVector(my_graph.nb_vertices, my_graph.vsize.constView(), my_buffer.view());
213
214 MPI_Gatherv(my_buffer.data(), my_graph.nb_vertices, MPI_INT,
215 buffer.data(), nb_vertices_per_rank.data(),
216 offset.data(), MPI_INT, io_rank, comm);
217
218 _convertVector(graph.nb_vertices, buffer.constView(), graph.vsize.view());
219 }
220
221 // Récupère la liste des numéros de sous-domaine d'un precedent partitionnement
222
223 if (is_master_io) {
224 graph.part.resize(graph.nb_vertices);
225 }
226
227 if (need_part) {
228
229 _convertVector(my_graph.nb_vertices, my_graph.part.constView(), my_buffer.view());
230
231 MPI_Gatherv(my_buffer.data(), my_graph.nb_vertices, MPI_INT,
232 buffer.data(), nb_vertices_per_rank.data(),
233 offset.data(), MPI_INT, io_rank, comm);
234
235 _convertVector(graph.nb_vertices, buffer.constView(), graph.part.view());
236 }
237
238 // Récupère la liste des poids aux sommets du graph. Il peut y en avoir plusieurs (ncon).
239
240 if (is_master_io) {
241 graph.vwgt.resize(graph.nb_vertices * ncon);
242 for (auto& x : offset) {
243 x *= ncon;
244 }
245 }
246
247 _convertVector(my_graph.nb_vertices * ncon, my_graph.vwgt.constView(), my_buffer.view());
248
249 MPI_Gatherv(my_buffer.data(), my_graph.nb_vertices * ncon, MPI_INT,
250 buffer.data(), nb_vertice_cons_per_rank.data(),
251 offset.data(), MPI_INT, io_rank, comm);
252
253 _convertVector(graph.nb_vertices * ncon, buffer.constView(), graph.vwgt.view());
254}
255
256/*---------------------------------------------------------------------------*/
257/*---------------------------------------------------------------------------*/
258
259void MetisGraphGather::
261{
262 MPI_Comm comm = _getMPICommunicator(m_parallel_mng);
263
264 const Int32 nb_rank = m_nb_rank;
265 const bool is_master_io = m_parallel_mng->isMasterIO();
266 int io_rank = 0;
267
268 UniqueArray<Int32> nb_vertices(nb_rank);
269 UniqueArray<Int32> displ(nb_rank);
270
271 for (int rank = 0; rank < nb_rank; ++rank) {
272 displ[rank] = CheckedConvert::toInt32(vtxdist[rank]);
273 nb_vertices[rank] = CheckedConvert::toInt32(vtxdist[rank + 1] - vtxdist[rank]);
274 }
275
278
279 if (is_master_io) {
280 send_buffer.resize(part.size());
281 }
282
283 _convertVector(send_buffer.size(), part, send_buffer.view());
284
285 MPI_Scatterv(send_buffer.data(), nb_vertices.data(),
286 displ.data(), MPI_INT, recv_buffer.data(),
287 my_part.size(), MPI_INT, io_rank, comm);
288
289 _convertVector(recv_buffer.size(), recv_buffer.constView(), my_part);
290}
291
292/*---------------------------------------------------------------------------*/
293/*---------------------------------------------------------------------------*/
294
295} // End namespace Arcane
296
297/*---------------------------------------------------------------------------*/
298/*---------------------------------------------------------------------------*/
#define ARCANE_FATAL(...)
Macro envoyant une exception FatalErrorException.
Interface du gestionnaire de parallélisme pour un sous-domaine.
virtual Parallel::Communicator communicator() const =0
Communicateur MPI associé à ce gestionnaire.
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Definition Lima.cc:120
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-