Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Don't use this anymore, as we are in the same namespace already.
[simgrid.git] / src / smpi / smpi_group.cpp
1 /* Copyright (c) 2010, 2013-2015. The SimGrid Team.
2  * All rights reserved.                                                     */
3
4 /* This program is free software; you can redistribute it and/or modify it
5  * under the terms of the license (GNU LGPL) which comes with this package. */
6
7 #include "private.h"
8
9 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_group, smpi, "Logging specific to SMPI (group)");
10  Group mpi_MPI_GROUP_EMPTY;
11 MPI_Group MPI_GROUP_EMPTY=&mpi_MPI_GROUP_EMPTY;
12
13 namespace simgrid{
14 namespace smpi{
15
16 Group::Group()
17 {
18   size_=0;                            /* size */
19   rank_to_index_map_=nullptr;                         /* rank_to_index_map_ */
20   index_to_rank_map_=nullptr;                         /* index_to_rank_map_ */
21   refcount_=1;                            /* refcount_: start > 0 so that this group never gets freed */
22 }
23
24 Group::Group(int n) : size_(n)
25 {
26   int i;
27   rank_to_index_map_ = xbt_new(int, size_);
28   index_to_rank_map_ = xbt_dict_new_homogeneous(xbt_free_f);
29   refcount_ = 1;
30   for (i = 0; i < size_; i++) {
31     rank_to_index_map_[i] = MPI_UNDEFINED;
32   }
33 }
34
35 Group::Group(MPI_Group origin)
36 {
37   char *key;
38   char *ptr_rank;
39   xbt_dict_cursor_t cursor = nullptr;
40   
41   int i;
42   if(origin != MPI_GROUP_NULL
43             && origin != MPI_GROUP_EMPTY)
44     {
45       size_ = origin->size();
46       rank_to_index_map_ = xbt_new(int, size_);
47       index_to_rank_map_ = xbt_dict_new_homogeneous(xbt_free_f);
48       refcount_ = 1;
49       for (i = 0; i < size_; i++) {
50         rank_to_index_map_[i] = origin->rank_to_index_map_[i];
51       }
52
53       xbt_dict_foreach(origin->index_to_rank_map_, cursor, key, ptr_rank) {
54         int * cp = static_cast<int*>(xbt_malloc(sizeof(int)));
55         *cp=*reinterpret_cast<int*>(ptr_rank);
56         xbt_dict_set(index_to_rank_map_, key, cp, nullptr);
57       }
58     }
59 }
60
61 Group::~Group()
62 {
63   xbt_free(rank_to_index_map_);
64   xbt_dict_free(&index_to_rank_map_);
65 }
66
67 void Group::destroy()
68 {
69   if(this != MPI_COMM_WORLD->group()
70           && this != MPI_GROUP_EMPTY)
71   this->unuse();
72 }
73
74 void Group::set_mapping(int index, int rank)
75 {
76   int * val_rank;
77
78   if (rank < size_) {
79     rank_to_index_map_[rank] = index;
80     if (index!=MPI_UNDEFINED ) {
81       val_rank = static_cast<int *>(xbt_malloc(sizeof(int)));
82       *val_rank = rank;
83
84       char * key = bprintf("%d", index);
85       xbt_dict_set(index_to_rank_map_, key, val_rank, nullptr);
86       xbt_free(key);
87     }
88   }
89 }
90
91 int Group::index(int rank)
92 {
93   int index = MPI_UNDEFINED;
94
95   if (0 <= rank && rank < size_) {
96     index = rank_to_index_map_[rank];
97   }
98   return index;
99 }
100
101 int Group::rank(int index)
102 {
103   int * ptr_rank = nullptr;
104   if (this==MPI_GROUP_EMPTY)
105     return MPI_UNDEFINED;
106   char * key = bprintf("%d", index);
107   ptr_rank = static_cast<int*>(xbt_dict_get_or_null(index_to_rank_map_, key));
108   xbt_free(key);
109
110   if (ptr_rank==nullptr)
111     return MPI_UNDEFINED;
112   return *ptr_rank;
113 }
114
115 int Group::use()
116 {
117   refcount_++;
118   return refcount_;
119 }
120
121 int Group::unuse()
122 {
123   refcount_--;
124   if (refcount_ <= 0) {
125     delete this;
126     return 0;
127   }
128   return refcount_;
129 }
130
131 int Group::size()
132 {
133   return size_;
134 }
135
136 int Group::compare(MPI_Group group2)
137 {
138   int result;
139   int i;
140   int index;
141   int rank;
142   int sz;
143
144   result = MPI_IDENT;
145   if (size_ != group2->size()) {
146     result = MPI_UNEQUAL;
147   } else {
148     sz = group2->size();
149     for (i = 0; i < sz; i++) {
150       index = this->index(i);
151       rank = group2->rank(index);
152       if (rank == MPI_UNDEFINED) {
153         result = MPI_UNEQUAL;
154         break;
155       }
156       if (rank != i) {
157         result = MPI_SIMILAR;
158       }
159     }
160   }
161   return result;
162 }
163
164 int Group::incl(int n, int* ranks, MPI_Group* newgroup)
165 {
166   int i=0;
167   int index=0;
168   if (n == 0) {
169     *newgroup = MPI_GROUP_EMPTY;
170   } else if (n == size_) {
171     *newgroup = this;
172     if(this!= MPI_COMM_WORLD->group()
173               && this != MPI_COMM_SELF->group()
174               && this != MPI_GROUP_EMPTY)
175     this->use();
176   } else {
177     *newgroup = new Group(n);
178     for (i = 0; i < n; i++) {
179       index = this->index(ranks[i]);
180       (*newgroup)->set_mapping(index, i);
181     }
182   }
183   return MPI_SUCCESS;
184 }
185
186 int Group::group_union(MPI_Group group2, MPI_Group* newgroup)
187 {
188   int size1 = size_;
189   int size2 = group2->size();
190   for (int i = 0; i < size2; i++) {
191     int proc2 = group2->index(i);
192     int proc1 = this->rank(proc2);
193     if (proc1 == MPI_UNDEFINED) {
194       size1++;
195     }
196   }
197   if (size1 == 0) {
198     *newgroup = MPI_GROUP_EMPTY;
199   } else {
200     *newgroup = new  Group(size1);
201     size2 = this->size();
202     for (int i = 0; i < size2; i++) {
203       int proc1 = this->index(i);
204       (*newgroup)->set_mapping(proc1, i);
205     }
206     for (int i = size2; i < size1; i++) {
207       int proc2 = group2->index(i - size2);
208       (*newgroup)->set_mapping(proc2, i);
209     }
210   }
211   return MPI_SUCCESS;
212 }
213
214 int Group::intersection(MPI_Group group2, MPI_Group* newgroup)
215 {
216   int size2 = group2->size();
217   for (int i = 0; i < size2; i++) {
218     int proc2 = group2->index(i);
219     int proc1 = this->rank(proc2);
220     if (proc1 == MPI_UNDEFINED) {
221       size2--;
222     }
223   }
224   if (size2 == 0) {
225     *newgroup = MPI_GROUP_EMPTY;
226   } else {
227     *newgroup = new  Group(size2);
228     int j=0;
229     for (int i = 0; i < group2->size(); i++) {
230       int proc2 = group2->index(i);
231       int proc1 = this->rank(proc2);
232       if (proc1 != MPI_UNDEFINED) {
233         (*newgroup)->set_mapping(proc2, j);
234         j++;
235       }
236     }
237   }
238   return MPI_SUCCESS;
239 }
240
241 int Group::difference(MPI_Group group2, MPI_Group* newgroup)
242 {
243   int newsize = size_;
244   int size2 = size_;
245   for (int i = 0; i < size2; i++) {
246     int proc1 = this->index(i);
247     int proc2 = group2->rank(proc1);
248     if (proc2 != MPI_UNDEFINED) {
249       newsize--;
250     }
251   }
252   if (newsize == 0) {
253     *newgroup = MPI_GROUP_EMPTY;
254   } else {
255     *newgroup = new  Group(newsize);
256     for (int i = 0; i < size2; i++) {
257       int proc1 = this->index(i);
258       int proc2 = group2->rank(proc1);
259       if (proc2 == MPI_UNDEFINED) {
260         (*newgroup)->set_mapping(proc1, i);
261       }
262     }
263   }
264   return MPI_SUCCESS;
265 }
266
267 int Group::excl(int n, int *ranks, MPI_Group * newgroup){
268   int oldsize = size_;
269   int newsize = oldsize - n;
270   *newgroup = new  Group(newsize);
271   int* to_exclude=xbt_new0(int, size_);
272   for (int i     = 0; i < oldsize; i++)
273     to_exclude[i]=0;
274   for (int i            = 0; i < n; i++)
275     to_exclude[ranks[i]]=1;
276   int j = 0;
277   for (int i = 0; i < oldsize; i++) {
278     if(to_exclude[i]==0){
279       int index = this->index(i);
280       (*newgroup)->set_mapping(index, j);
281       j++;
282     }
283   }
284   xbt_free(to_exclude);
285   return MPI_SUCCESS;
286
287 }
288
289 int Group::range_incl(int n, int ranges[][3], MPI_Group * newgroup){
290   int newsize = 0;
291   for (int i = 0; i < n; i++) {
292     for (int rank = ranges[i][0];                    /* First */
293          rank >= 0 && rank < size_; /* Last */
294          ) {
295       newsize++;
296       if(rank == ranges[i][1]){/*already last ?*/
297         break;
298       }
299       rank += ranges[i][2]; /* Stride */
300       if (ranges[i][0] < ranges[i][1]) {
301         if (rank > ranges[i][1])
302           break;
303       } else {
304         if (rank < ranges[i][1])
305           break;
306       }
307     }
308   }
309   *newgroup = new  Group(newsize);
310   int j     = 0;
311   for (int i = 0; i < n; i++) {
312     for (int rank = ranges[i][0];                    /* First */
313          rank >= 0 && rank < size_; /* Last */
314          ) {
315       int index = this->index(rank);
316       (*newgroup)->set_mapping(index, j);
317       j++;
318       if(rank == ranges[i][1]){/*already last ?*/
319         break;
320       }
321       rank += ranges[i][2]; /* Stride */
322       if (ranges[i][0] < ranges[i][1]) {
323         if (rank > ranges[i][1])
324           break;
325       } else {
326         if (rank < ranges[i][1])
327           break;
328       }
329     }
330   }
331   return MPI_SUCCESS;
332 }
333
334 int Group::range_excl(int n, int ranges[][3], MPI_Group * newgroup){
335   int newsize = size_;
336   for (int i = 0; i < n; i++) {
337     for (int rank = ranges[i][0];                    /* First */
338          rank >= 0 && rank < size_; /* Last */
339          ) {
340       newsize--;
341       if(rank == ranges[i][1]){/*already last ?*/
342         break;
343       }
344       rank += ranges[i][2]; /* Stride */
345       if (ranges[i][0] < ranges[i][1]) {
346         if (rank > ranges[i][1])
347           break;
348       } else {
349         if (rank < ranges[i][1])
350           break;
351       }
352     }
353   }
354   if (newsize == 0) {
355     *newgroup = MPI_GROUP_EMPTY;
356   } else {
357     *newgroup = new  Group(newsize);
358     int newrank = 0;
359     int oldrank = 0;
360     while (newrank < newsize) {
361       int add = 1;
362       for (int i = 0; i < n; i++) {
363         for (int rank = ranges[i][0]; rank >= 0 && rank < size_;) {
364           if(rank==oldrank){
365             add = 0;
366             break;
367           }
368           if(rank == ranges[i][1]){/*already last ?*/
369             break;
370           }
371           rank += ranges[i][2]; /* Stride */
372           if (ranges[i][0]<ranges[i][1]){
373             if (rank > ranges[i][1])
374               break;
375           }else{
376             if (rank < ranges[i][1])
377               break;
378           }
379         }
380       }
381       if(add==1){
382         int index = this->index(oldrank);
383         (*newgroup)->set_mapping(index, newrank);
384         newrank++;
385       }
386       oldrank++;
387     }
388   }
389   return MPI_SUCCESS;
390 }
391
392 }
393 }