Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
leak--
[simgrid.git] / src / smpi / colls / smpi_openmpi_selector.cpp
1 /* selector for collective algorithms based on openmpi's default coll_tuned_decision_fixed selector
2  * Updated 02/2022                                                          */
3
4 /* Copyright (c) 2009-2022. The SimGrid Team.
5  * All rights reserved.                                                     */
6
7 /* This program is free software; you can redistribute it and/or modify it
8  * under the terms of the license (GNU LGPL) which comes with this package. */
9
10 #include "colls_private.hpp"
11
12 #include <memory>
13
14 namespace simgrid {
15 namespace smpi {
16
17 int allreduce__ompi(const void *sbuf, void *rbuf, int count,
18                     MPI_Datatype dtype, MPI_Op op, MPI_Comm comm)
19 {
20     size_t total_dsize = dtype->size() * (ptrdiff_t)count;
21     int communicator_size = comm->size();
22     int alg = 1;
23     int(*funcs[]) (const void*, void*, int, MPI_Datatype, MPI_Op, MPI_Comm)={
24         &allreduce__redbcast,
25         &allreduce__redbcast,
26         &allreduce__rdb,
27         &allreduce__lr,
28         &allreduce__ompi_ring_segmented,
29         &allreduce__rab_rdb
30     };
31
32     /** Algorithms:
33      *  {1, "basic_linear"},
34      *  {2, "nonoverlapping"},
35      *  {3, "recursive_doubling"},
36      *  {4, "ring"},
37      *  {5, "segmented_ring"},
38      *  {6, "rabenseifner"
39      *
40      * Currently, ring, segmented ring, and rabenseifner do not support
41      * non-commutative operations.
42      */
43     if ((op != MPI_OP_NULL) && not op->is_commutative()) {
44         if (communicator_size < 4) {
45             if (total_dsize < 131072) {
46                 alg = 3;
47             } else {
48                 alg = 1;
49             }
50         } else if (communicator_size < 8) {
51             alg = 3;
52         } else if (communicator_size < 16) {
53             if (total_dsize < 1048576) {
54                 alg = 3;
55             } else {
56                 alg = 2;
57             }
58         } else if (communicator_size < 128) {
59             alg = 3;
60         } else if (communicator_size < 256) {
61             if (total_dsize < 131072) {
62                 alg = 2;
63             } else if (total_dsize < 524288) {
64                 alg = 3;
65             } else {
66                 alg = 2;
67             }
68         } else if (communicator_size < 512) {
69             if (total_dsize < 4096) {
70                 alg = 2;
71             } else if (total_dsize < 524288) {
72                 alg = 3;
73             } else {
74                 alg = 2;
75             }
76         } else {
77             if (total_dsize < 2048) {
78                 alg = 2;
79             } else {
80                 alg = 3;
81             }
82         }
83     } else {
84         if (communicator_size < 4) {
85             if (total_dsize < 8) {
86                 alg = 4;
87             } else if (total_dsize < 4096) {
88                 alg = 3;
89             } else if (total_dsize < 8192) {
90                 alg = 4;
91             } else if (total_dsize < 16384) {
92                 alg = 3;
93             } else if (total_dsize < 65536) {
94                 alg = 4;
95             } else if (total_dsize < 262144) {
96                 alg = 5;
97             } else {
98                 alg = 6;
99             }
100         } else if (communicator_size < 8) {
101             if (total_dsize < 16) {
102                 alg = 4;
103             } else if (total_dsize < 8192) {
104                 alg = 3;
105             } else {
106                 alg = 6;
107             }
108         } else if (communicator_size < 16) {
109             if (total_dsize < 8192) {
110                 alg = 3;
111             } else {
112                 alg = 6;
113             }
114         } else if (communicator_size < 32) {
115             if (total_dsize < 64) {
116                 alg = 5;
117             } else if (total_dsize < 4096) {
118                 alg = 3;
119             } else {
120                 alg = 6;
121             }
122         } else if (communicator_size < 64) {
123             if (total_dsize < 128) {
124                 alg = 5;
125             } else {
126                 alg = 6;
127             }
128         } else if (communicator_size < 128) {
129             if (total_dsize < 262144) {
130                 alg = 3;
131             } else {
132                 alg = 6;
133             }
134         } else if (communicator_size < 256) {
135             if (total_dsize < 131072) {
136                 alg = 2;
137             } else if (total_dsize < 262144) {
138                 alg = 3;
139             } else {
140                 alg = 6;
141             }
142         } else if (communicator_size < 512) {
143             if (total_dsize < 4096) {
144                 alg = 2;
145             } else {
146                 alg = 6;
147             }
148         } else if (communicator_size < 2048) {
149             if (total_dsize < 2048) {
150                 alg = 2;
151             } else if (total_dsize < 16384) {
152                 alg = 3;
153             } else {
154                 alg = 6;
155             }
156         } else if (communicator_size < 4096) {
157             if (total_dsize < 2048) {
158                 alg = 2;
159             } else if (total_dsize < 4096) {
160                 alg = 5;
161             } else if (total_dsize < 16384) {
162                 alg = 3;
163             } else {
164                 alg = 6;
165             }
166         } else {
167             if (total_dsize < 2048) {
168                 alg = 2;
169             } else if (total_dsize < 16384) {
170                 alg = 5;
171             } else if (total_dsize < 32768) {
172                 alg = 3;
173             } else {
174                 alg = 6;
175             }
176         }
177     }
178     return funcs[alg-1](sbuf, rbuf, count, dtype, op, comm);
179 }
180
181
182
183 int alltoall__ompi(const void *sbuf, int scount,
184                    MPI_Datatype sdtype,
185                    void* rbuf, int rcount,
186                    MPI_Datatype rdtype,
187                    MPI_Comm comm)
188 {
189     int alg = 1;
190     size_t dsize, total_dsize;
191     int communicator_size = comm->size();
192
193     if (MPI_IN_PLACE != sbuf) {
194         dsize = sdtype->size();
195     } else {
196         dsize = rdtype->size();
197     }
198     total_dsize = dsize * (ptrdiff_t)scount;
199     int (*funcs[])(const void *, int, MPI_Datatype, void*, int, MPI_Datatype, MPI_Comm) = {
200         &alltoall__basic_linear,
201         &alltoall__pair,
202         &alltoall__bruck,
203         &alltoall__basic_linear,
204         &alltoall__basic_linear
205     };
206     /** Algorithms:
207      *  {1, "linear"},
208      *  {2, "pairwise"},
209      *  {3, "modified_bruck"},
210      *  {4, "linear_sync"},
211      *  {5, "two_proc"},
212      */
213     if (communicator_size == 2) {
214         if (total_dsize < 2) {
215             alg = 2;
216         } else if (total_dsize < 4) {
217             alg = 5;
218         } else if (total_dsize < 16) {
219             alg = 2;
220         } else if (total_dsize < 64) {
221             alg = 5;
222         } else if (total_dsize < 256) {
223             alg = 2;
224         } else if (total_dsize < 4096) {
225             alg = 5;
226         } else if (total_dsize < 32768) {
227             alg = 2;
228         } else if (total_dsize < 262144) {
229             alg = 4;
230         } else if (total_dsize < 1048576) {
231             alg = 5;
232         } else {
233             alg = 2;
234         }
235     } else if (communicator_size < 8) {
236         if (total_dsize < 8192) {
237             alg = 4;
238         } else if (total_dsize < 16384) {
239             alg = 1;
240         } else if (total_dsize < 65536) {
241             alg = 4;
242         } else if (total_dsize < 524288) {
243             alg = 1;
244         } else if (total_dsize < 1048576) {
245             alg = 2;
246         } else {
247             alg = 1;
248         }
249     } else if (communicator_size < 16) {
250         if (total_dsize < 262144) {
251             alg = 4;
252         } else {
253             alg = 1;
254         }
255     } else if (communicator_size < 32) {
256         if (total_dsize < 4) {
257             alg = 4;
258         } else if (total_dsize < 512) {
259             alg = 3;
260         } else if (total_dsize < 8192) {
261             alg = 4;
262         } else if (total_dsize < 32768) {
263             alg = 1;
264         } else if (total_dsize < 262144) {
265             alg = 4;
266         } else if (total_dsize < 524288) {
267             alg = 1;
268         } else {
269             alg = 4;
270         }
271     } else if (communicator_size < 64) {
272         if (total_dsize < 512) {
273             alg = 3;
274         } else if (total_dsize < 524288) {
275             alg = 1;
276         } else {
277             alg = 4;
278         }
279     } else if (communicator_size < 128) {
280         if (total_dsize < 1024) {
281             alg = 3;
282         } else if (total_dsize < 2048) {
283             alg = 1;
284         } else if (total_dsize < 4096) {
285             alg = 4;
286         } else if (total_dsize < 262144) {
287             alg = 1;
288         } else {
289             alg = 2;
290         }
291     } else if (communicator_size < 256) {
292         if (total_dsize < 1024) {
293             alg = 3;
294         } else if (total_dsize < 2048) {
295             alg = 4;
296         } else if (total_dsize < 262144) {
297             alg = 1;
298         } else {
299             alg = 2;
300         }
301     } else if (communicator_size < 512) {
302         if (total_dsize < 1024) {
303             alg = 3;
304         } else if (total_dsize < 8192) {
305             alg = 4;
306         } else if (total_dsize < 32768) {
307             alg = 1;
308         } else {
309             alg = 2;
310         }
311     } else if (communicator_size < 1024) {
312         if (total_dsize < 512) {
313             alg = 3;
314         } else if (total_dsize < 8192) {
315             alg = 4;
316         } else if (total_dsize < 16384) {
317             alg = 1;
318         } else if (total_dsize < 131072) {
319             alg = 4;
320         } else if (total_dsize < 262144) {
321             alg = 1;
322         } else {
323             alg = 2;
324         }
325     } else if (communicator_size < 2048) {
326         if (total_dsize < 512) {
327             alg = 3;
328         } else if (total_dsize < 1024) {
329             alg = 4;
330         } else if (total_dsize < 2048) {
331             alg = 1;
332         } else if (total_dsize < 16384) {
333             alg = 4;
334         } else if (total_dsize < 262144) {
335             alg = 1;
336         } else {
337             alg = 4;
338         }
339     } else if (communicator_size < 4096) {
340         if (total_dsize < 1024) {
341             alg = 3;
342         } else if (total_dsize < 4096) {
343             alg = 4;
344         } else if (total_dsize < 8192) {
345             alg = 1;
346         } else if (total_dsize < 131072) {
347             alg = 4;
348         } else {
349             alg = 1;
350         }
351     } else {
352         if (total_dsize < 2048) {
353             alg = 3;
354         } else if (total_dsize < 8192) {
355             alg = 4;
356         } else if (total_dsize < 16384) {
357             alg = 1;
358         } else if (total_dsize < 32768) {
359             alg = 4;
360         } else if (total_dsize < 65536) {
361             alg = 1;
362         } else {
363             alg = 4;
364         }
365     }
366
367     return funcs[alg-1](sbuf, scount, sdtype,
368                           rbuf, rcount, rdtype, comm);
369 }
370
371 int alltoallv__ompi(const void *sbuf, const int *scounts, const int *sdisps,
372                     MPI_Datatype sdtype,
373                     void *rbuf, const int *rcounts, const int *rdisps,
374                     MPI_Datatype rdtype,
375                     MPI_Comm  comm
376                     )
377 {
378     int communicator_size = comm->size();
379     int alg = 1;
380     int (*funcs[])(const void *, const int*, const int*, MPI_Datatype, void*, const int*, const int*, MPI_Datatype, MPI_Comm) = {
381         &alltoallv__ompi_basic_linear,
382         &alltoallv__pair
383     };
384    /** Algorithms:
385      *  {1, "basic_linear"},
386      *  {2, "pairwise"},
387      *
388      * We can only optimize based on com size
389      */
390     if (communicator_size < 4) {
391         alg = 2;
392     } else if (communicator_size < 64) {
393         alg = 1;
394     } else if (communicator_size < 128) {
395         alg = 2;
396     } else if (communicator_size < 256) {
397         alg = 1;
398     } else if (communicator_size < 1024) {
399         alg = 2;
400     } else {
401         alg = 1;
402     }
403     return funcs[alg-1](sbuf, scounts, sdisps, sdtype,
404                            rbuf, rcounts, rdisps,rdtype,
405                            comm);
406 }
407
408 int barrier__ompi(MPI_Comm  comm)
409 {
410     int communicator_size = comm->size();
411     int alg = 1;
412     int (*funcs[])(MPI_Comm) = {
413         &barrier__ompi_basic_linear,
414         &barrier__ompi_basic_linear,
415         &barrier__ompi_recursivedoubling,
416         &barrier__ompi_bruck,
417         &barrier__ompi_two_procs,
418         &barrier__ompi_tree
419     };
420     /** Algorithms:
421      *  {1, "linear"},
422      *  {2, "double_ring"},
423      *  {3, "recursive_doubling"},
424      *  {4, "bruck"},
425      *  {5, "two_proc"},
426      *  {6, "tree"},
427      *
428      * We can only optimize based on com size
429      */
430     if (communicator_size < 4) {
431         alg = 3;
432     } else if (communicator_size < 8) {
433         alg = 1;
434     } else if (communicator_size < 64) {
435         alg = 3;
436     } else if (communicator_size < 256) {
437         alg = 4;
438     } else if (communicator_size < 512) {
439         alg = 6;
440     } else if (communicator_size < 1024) {
441         alg = 4;
442     } else if (communicator_size < 4096) {
443         alg = 6;
444     } else {
445         alg = 4;
446     }
447
448     return funcs[alg-1](comm);
449 }
450
451 int bcast__ompi(void *buff, int count, MPI_Datatype datatype, int root, MPI_Comm  comm)
452 {
453     int alg = 1;
454     size_t total_dsize, dsize;
455
456     int communicator_size = comm->size();
457
458     dsize = datatype->size();
459     total_dsize = dsize * (unsigned long)count;
460     int (*funcs[])(void*, int, MPI_Datatype, int, MPI_Comm) = {
461         &bcast__NTSL,
462         &bcast__ompi_pipeline,
463         &bcast__ompi_pipeline,
464         &bcast__ompi_split_bintree,
465         &bcast__NTSB,
466         &bcast__binomial_tree,
467         &bcast__mvapich2_knomial_intra_node,
468         &bcast__scatter_rdb_allgather,
469         &bcast__scatter_LR_allgather,
470     };
471     /** Algorithms:
472      *  {1, "basic_linear"},
473      *  {2, "chain"},
474      *  {3, "pipeline"},
475      *  {4, "split_binary_tree"},
476      *  {5, "binary_tree"},
477      *  {6, "binomial"},
478      *  {7, "knomial"},
479      *  {8, "scatter_allgather"},
480      *  {9, "scatter_allgather_ring"},
481      */
482     if (communicator_size < 4) {
483         if (total_dsize < 32) {
484             alg = 3;
485         } else if (total_dsize < 256) {
486             alg = 5;
487         } else if (total_dsize < 512) {
488             alg = 3;
489         } else if (total_dsize < 1024) {
490             alg = 7;
491         } else if (total_dsize < 32768) {
492             alg = 1;
493         } else if (total_dsize < 131072) {
494             alg = 5;
495         } else if (total_dsize < 262144) {
496             alg = 2;
497         } else if (total_dsize < 524288) {
498             alg = 1;
499         } else if (total_dsize < 1048576) {
500             alg = 6;
501         } else {
502             alg = 5;
503         }
504     } else if (communicator_size < 8) {
505         if (total_dsize < 64) {
506             alg = 5;
507         } else if (total_dsize < 128) {
508             alg = 6;
509         } else if (total_dsize < 2048) {
510             alg = 5;
511         } else if (total_dsize < 8192) {
512             alg = 6;
513         } else if (total_dsize < 1048576) {
514             alg = 1;
515         } else {
516             alg = 2;
517         }
518     } else if (communicator_size < 16) {
519         if (total_dsize < 8) {
520             alg = 7;
521         } else if (total_dsize < 64) {
522             alg = 5;
523         } else if (total_dsize < 4096) {
524             alg = 7;
525         } else if (total_dsize < 16384) {
526             alg = 5;
527         } else if (total_dsize < 32768) {
528             alg = 6;
529         } else {
530             alg = 1;
531         }
532     } else if (communicator_size < 32) {
533         if (total_dsize < 4096) {
534             alg = 7;
535         } else if (total_dsize < 1048576) {
536             alg = 6;
537         } else {
538             alg = 8;
539         }
540     } else if (communicator_size < 64) {
541         if (total_dsize < 2048) {
542             alg = 6;
543         } else {
544             alg = 7;
545         }
546     } else if (communicator_size < 128) {
547         alg = 7;
548     } else if (communicator_size < 256) {
549         if (total_dsize < 2) {
550             alg = 6;
551         } else if (total_dsize < 16384) {
552             alg = 5;
553         } else if (total_dsize < 32768) {
554             alg = 1;
555         } else if (total_dsize < 65536) {
556             alg = 5;
557         } else {
558             alg = 7;
559         }
560     } else if (communicator_size < 1024) {
561         if (total_dsize < 16384) {
562             alg = 7;
563         } else if (total_dsize < 32768) {
564             alg = 4;
565         } else {
566             alg = 7;
567         }
568     } else if (communicator_size < 2048) {
569         if (total_dsize < 524288) {
570             alg = 7;
571         } else {
572             alg = 8;
573         }
574     } else if (communicator_size < 4096) {
575         if (total_dsize < 262144) {
576             alg = 7;
577         } else {
578             alg = 8;
579         }
580     } else {
581         if (total_dsize < 8192) {
582             alg = 7;
583         } else if (total_dsize < 16384) {
584             alg = 5;
585         } else if (total_dsize < 262144) {
586             alg = 7;
587         } else {
588             alg = 8;
589         }
590     }
591     return funcs[alg-1](buff, count, datatype, root, comm);
592 }
593
594 int reduce__ompi(const void *sendbuf, void *recvbuf,
595                  int count, MPI_Datatype  datatype,
596                  MPI_Op   op, int root,
597                  MPI_Comm   comm)
598 {
599     size_t total_dsize, dsize;
600     int alg = 1;
601     int communicator_size = comm->size();
602
603     dsize=datatype->size();
604     total_dsize = dsize * count;
605     int (*funcs[])(const void*, void*, int, MPI_Datatype, MPI_Op, int, MPI_Comm) = {
606         &reduce__ompi_basic_linear,
607         &reduce__ompi_chain,
608         &reduce__ompi_pipeline,
609         &reduce__ompi_binary,
610         &reduce__ompi_binomial,
611         &reduce__ompi_in_order_binary,
612         //&reduce__rab our rab can't be used with all datatypes
613         &reduce__ompi_basic_linear
614     };
615     /** Algorithms:
616      *  {1, "linear"},
617      *  {2, "chain"},
618      *  {3, "pipeline"},
619      *  {4, "binary"},
620      *  {5, "binomial"},
621      *  {6, "in-order_binary"},
622      *  {7, "rabenseifner"},
623      *
624      * Currently, only linear and in-order binary tree algorithms are
625      * capable of non commutative ops.
626      */
627      if ((op != MPI_OP_NULL) && not op->is_commutative()) {
628         if (communicator_size < 4) {
629             if (total_dsize < 8) {
630                 alg = 6;
631             } else {
632                 alg = 1;
633             }
634         } else if (communicator_size < 8) {
635             alg = 1;
636         } else if (communicator_size < 16) {
637             if (total_dsize < 1024) {
638                 alg = 6;
639             } else if (total_dsize < 8192) {
640                 alg = 1;
641             } else if (total_dsize < 16384) {
642                 alg = 6;
643             } else if (total_dsize < 262144) {
644                 alg = 1;
645             } else {
646                 alg = 6;
647             }
648         } else if (communicator_size < 128) {
649             alg = 6;
650         } else if (communicator_size < 256) {
651             if (total_dsize < 512) {
652                 alg = 6;
653             } else if (total_dsize < 1024) {
654                 alg = 1;
655             } else {
656                 alg = 6;
657             }
658         } else {
659             alg = 6;
660         }
661     } else {
662         if (communicator_size < 4) {
663             if (total_dsize < 8) {
664                 alg = 7;
665             } else if (total_dsize < 16) {
666                 alg = 4;
667             } else if (total_dsize < 32) {
668                 alg = 3;
669             } else if (total_dsize < 262144) {
670                 alg = 1;
671             } else if (total_dsize < 524288) {
672                 alg = 3;
673             } else if (total_dsize < 1048576) {
674                 alg = 2;
675             } else {
676                 alg = 3;
677             }
678         } else if (communicator_size < 8) {
679             if (total_dsize < 4096) {
680                 alg = 4;
681             } else if (total_dsize < 65536) {
682                 alg = 2;
683             } else if (total_dsize < 262144) {
684                 alg = 5;
685             } else if (total_dsize < 524288) {
686                 alg = 1;
687             } else if (total_dsize < 1048576) {
688                 alg = 5;
689             } else {
690                 alg = 1;
691             }
692         } else if (communicator_size < 16) {
693             if (total_dsize < 8192) {
694                 alg = 4;
695             } else {
696                 alg = 5;
697             }
698         } else if (communicator_size < 32) {
699             if (total_dsize < 4096) {
700                 alg = 4;
701             } else {
702                 alg = 5;
703             }
704         } else if (communicator_size < 256) {
705             alg = 5;
706         } else if (communicator_size < 512) {
707             if (total_dsize < 8192) {
708                 alg = 5;
709             } else if (total_dsize < 16384) {
710                 alg = 6;
711             } else {
712                 alg = 5;
713             }
714         } else if (communicator_size < 2048) {
715             alg = 5;
716         } else if (communicator_size < 4096) {
717             if (total_dsize < 512) {
718                 alg = 5;
719             } else if (total_dsize < 1024) {
720                 alg = 6;
721             } else if (total_dsize < 8192) {
722                 alg = 5;
723             } else if (total_dsize < 16384) {
724                 alg = 6;
725             } else {
726                 alg = 5;
727             }
728         } else {
729             if (total_dsize < 16) {
730                 alg = 5;
731             } else if (total_dsize < 32) {
732                 alg = 6;
733             } else if (total_dsize < 1024) {
734                 alg = 5;
735             } else if (total_dsize < 2048) {
736                 alg = 6;
737             } else if (total_dsize < 8192) {
738                 alg = 5;
739             } else if (total_dsize < 16384) {
740                 alg = 6;
741             } else {
742                 alg = 5;
743             }
744         }
745     }
746
747     return funcs[alg-1] (sendbuf, recvbuf, count, datatype, op, root, comm);
748 }
749
750 int reduce_scatter__ompi(const void *sbuf, void *rbuf,
751                          const int *rcounts,
752                          MPI_Datatype dtype,
753                          MPI_Op  op,
754                          MPI_Comm  comm
755                          )
756 {
757     size_t total_dsize, dsize;
758     int communicator_size = comm->size();
759     int alg = 1;
760     int zerocounts = 0;
761     dsize=dtype->size();
762     total_dsize = 0;
763     for (int i = 0; i < communicator_size; i++) {
764         total_dsize += rcounts[i];
765        // if (0 == rcounts[i]) {
766         //    zerocounts = 1;
767         //}
768     }
769     total_dsize *= dsize;
770     int (*funcs[])(const void*, void*, const int*, MPI_Datatype, MPI_Op, MPI_Comm) = {
771         &reduce_scatter__default,
772         &reduce_scatter__ompi_basic_recursivehalving,
773         &reduce_scatter__ompi_ring,
774         &reduce_scatter__ompi_butterfly,
775     };
776     /** Algorithms:
777      *  {1, "non-overlapping"},
778      *  {2, "recursive_halving"},
779      *  {3, "ring"},
780      *  {4, "butterfly"},
781      *
782      * Non commutative algorithm capability needs re-investigation.
783      * Defaulting to non overlapping for non commutative ops.
784      */
785     if (((op != MPI_OP_NULL) && not op->is_commutative()) || (zerocounts)) {
786         alg = 1;
787     } else {
788         if (communicator_size < 4) {
789             if (total_dsize < 65536) {
790                 alg = 3;
791             } else if (total_dsize < 131072) {
792                 alg = 4;
793             } else {
794                 alg = 3;
795             }
796         } else if (communicator_size < 8) {
797             if (total_dsize < 8) {
798                 alg = 1;
799             } else if (total_dsize < 262144) {
800                 alg = 2;
801             } else {
802                 alg = 3;
803             }
804         } else if (communicator_size < 32) {
805             if (total_dsize < 262144) {
806                 alg = 2;
807             } else {
808                 alg = 3;
809             }
810         } else if (communicator_size < 64) {
811             if (total_dsize < 64) {
812                 alg = 1;
813             } else if (total_dsize < 2048) {
814                 alg = 2;
815             } else if (total_dsize < 524288) {
816                 alg = 4;
817             } else {
818                 alg = 3;
819             }
820         } else if (communicator_size < 128) {
821             if (total_dsize < 256) {
822                 alg = 1;
823             } else if (total_dsize < 512) {
824                 alg = 2;
825             } else if (total_dsize < 2048) {
826                 alg = 4;
827             } else if (total_dsize < 4096) {
828                 alg = 2;
829             } else {
830                 alg = 4;
831             }
832         } else if (communicator_size < 256) {
833             if (total_dsize < 256) {
834                 alg = 1;
835             } else if (total_dsize < 512) {
836                 alg = 2;
837             } else {
838                 alg = 4;
839             }
840         } else if (communicator_size < 512) {
841             if (total_dsize < 256) {
842                 alg = 1;
843             } else if (total_dsize < 1024) {
844                 alg = 2;
845             } else {
846                 alg = 4;
847             }
848         } else if (communicator_size < 1024) {
849             if (total_dsize < 512) {
850                 alg = 1;
851             } else if (total_dsize < 2048) {
852                 alg = 2;
853             } else if (total_dsize < 8192) {
854                 alg = 4;
855             } else if (total_dsize < 16384) {
856                 alg = 2;
857             } else {
858                 alg = 4;
859             }
860         } else if (communicator_size < 2048) {
861             if (total_dsize < 512) {
862                 alg = 1;
863             } else if (total_dsize < 4096) {
864                 alg = 2;
865             } else if (total_dsize < 16384) {
866                 alg = 4;
867             } else if (total_dsize < 32768) {
868                 alg = 2;
869             } else {
870                 alg = 4;
871             }
872         } else if (communicator_size < 4096) {
873             if (total_dsize < 512) {
874                 alg = 1;
875             } else if (total_dsize < 4096) {
876                 alg = 2;
877             } else {
878                 alg = 4;
879             }
880         } else {
881             if (total_dsize < 1024) {
882                 alg = 1;
883             } else if (total_dsize < 8192) {
884                 alg = 2;
885             } else {
886                 alg = 4;
887             }
888         }
889     }
890
891     return funcs[alg-1] (sbuf, rbuf, rcounts, dtype, op, comm);
892 }
893
894 int allgather__ompi(const void *sbuf, int scount,
895                     MPI_Datatype sdtype,
896                     void* rbuf, int rcount,
897                     MPI_Datatype rdtype,
898                     MPI_Comm  comm
899                     )
900 {
901     int communicator_size;
902     size_t dsize, total_dsize;
903     int alg = 1;
904     communicator_size = comm->size();
905     if (MPI_IN_PLACE != sbuf) {
906         dsize = sdtype->size();
907     } else {
908         dsize = rdtype->size();
909     }
910     total_dsize = dsize * (ptrdiff_t)scount;
911     int (*funcs[])(const void*, int, MPI_Datatype, void*, int, MPI_Datatype, MPI_Comm) = {
912         &allgather__NTSLR_NB,
913         &allgather__bruck,
914         &allgather__rdb,
915         &allgather__ring,
916         &allgather__ompi_neighborexchange,
917         &allgather__pair
918     };
919     /** Algorithms:
920      *  {1, "linear"},
921      *  {2, "bruck"},
922      *  {3, "recursive_doubling"},
923      *  {4, "ring"},
924      *  {5, "neighbor"},
925      *  {6, "two_proc"}
926      */
927     if (communicator_size == 2) {
928         alg = 6;
929     } else if (communicator_size < 32) {
930         alg = 3;
931     } else if (communicator_size < 64) {
932         if (total_dsize < 1024) {
933             alg = 3;
934         } else if (total_dsize < 65536) {
935             alg = 5;
936         } else {
937             alg = 4;
938         }
939     } else if (communicator_size < 128) {
940         if (total_dsize < 512) {
941             alg = 3;
942         } else if (total_dsize < 65536) {
943             alg = 5;
944         } else {
945             alg = 4;
946         }
947     } else if (communicator_size < 256) {
948         if (total_dsize < 512) {
949             alg = 3;
950         } else if (total_dsize < 131072) {
951             alg = 5;
952         } else if (total_dsize < 524288) {
953             alg = 4;
954         } else if (total_dsize < 1048576) {
955             alg = 5;
956         } else {
957             alg = 4;
958         }
959     } else if (communicator_size < 512) {
960         if (total_dsize < 32) {
961             alg = 3;
962         } else if (total_dsize < 128) {
963             alg = 2;
964         } else if (total_dsize < 1024) {
965             alg = 3;
966         } else if (total_dsize < 131072) {
967             alg = 5;
968         } else if (total_dsize < 524288) {
969             alg = 4;
970         } else if (total_dsize < 1048576) {
971             alg = 5;
972         } else {
973             alg = 4;
974         }
975     } else if (communicator_size < 1024) {
976         if (total_dsize < 64) {
977             alg = 3;
978         } else if (total_dsize < 256) {
979             alg = 2;
980         } else if (total_dsize < 2048) {
981             alg = 3;
982         } else {
983             alg = 5;
984         }
985     } else if (communicator_size < 2048) {
986         if (total_dsize < 4) {
987             alg = 3;
988         } else if (total_dsize < 8) {
989             alg = 2;
990         } else if (total_dsize < 16) {
991             alg = 3;
992         } else if (total_dsize < 32) {
993             alg = 2;
994         } else if (total_dsize < 256) {
995             alg = 3;
996         } else if (total_dsize < 512) {
997             alg = 2;
998         } else if (total_dsize < 4096) {
999             alg = 3;
1000         } else {
1001             alg = 5;
1002         }
1003     } else if (communicator_size < 4096) {
1004         if (total_dsize < 32) {
1005             alg = 2;
1006         } else if (total_dsize < 128) {
1007             alg = 3;
1008         } else if (total_dsize < 512) {
1009             alg = 2;
1010         } else if (total_dsize < 4096) {
1011             alg = 3;
1012         } else {
1013             alg = 5;
1014         }
1015     } else {
1016         if (total_dsize < 2) {
1017             alg = 3;
1018         } else if (total_dsize < 8) {
1019             alg = 2;
1020         } else if (total_dsize < 16) {
1021             alg = 3;
1022         } else if (total_dsize < 512) {
1023             alg = 2;
1024         } else if (total_dsize < 4096) {
1025             alg = 3;
1026         } else {
1027             alg = 5;
1028         }
1029     }
1030
1031     return funcs[alg-1](sbuf, scount, sdtype, rbuf, rcount, rdtype, comm);
1032
1033 }
1034
1035 int allgatherv__ompi(const void *sbuf, int scount,
1036                      MPI_Datatype sdtype,
1037                      void* rbuf, const int *rcounts,
1038                      const int *rdispls,
1039                      MPI_Datatype rdtype,
1040                      MPI_Comm  comm
1041                      )
1042 {
1043     int i;
1044     int communicator_size;
1045     size_t dsize, total_dsize;
1046     int alg = 1;
1047     communicator_size = comm->size();
1048     if (MPI_IN_PLACE != sbuf) {
1049         dsize = sdtype->size();
1050     } else {
1051         dsize = rdtype->size();
1052     }
1053
1054     total_dsize = 0;
1055     for (i = 0; i < communicator_size; i++) {
1056         total_dsize += dsize * rcounts[i];
1057     }
1058
1059     /* use the per-rank data size as basis, similar to allgather */
1060     size_t per_rank_dsize = total_dsize / communicator_size;
1061
1062     int (*funcs[])(const void*, int, MPI_Datatype, void*, const int*, const int*, MPI_Datatype, MPI_Comm) = {
1063         &allgatherv__GB,
1064         &allgatherv__ompi_bruck,
1065         &allgatherv__mpich_ring,
1066         &allgatherv__ompi_neighborexchange,
1067         &allgatherv__pair
1068     };
1069     /** Algorithms:
1070      *  {1, "default"},
1071      *  {2, "bruck"},
1072      *  {3, "ring"},
1073      *  {4, "neighbor"},
1074      *  {5, "two_proc"},
1075      */
1076     if (communicator_size == 2) {
1077         if (per_rank_dsize < 2048) {
1078             alg = 3;
1079         } else if (per_rank_dsize < 4096) {
1080             alg = 5;
1081         } else if (per_rank_dsize < 8192) {
1082             alg = 3;
1083         } else {
1084             alg = 5;
1085         }
1086     } else if (communicator_size < 8) {
1087         if (per_rank_dsize < 256) {
1088             alg = 1;
1089         } else if (per_rank_dsize < 4096) {
1090             alg = 4;
1091         } else if (per_rank_dsize < 8192) {
1092             alg = 3;
1093         } else if (per_rank_dsize < 16384) {
1094             alg = 4;
1095         } else if (per_rank_dsize < 262144) {
1096             alg = 2;
1097         } else {
1098             alg = 4;
1099         }
1100     } else if (communicator_size < 16) {
1101         if (per_rank_dsize < 1024) {
1102             alg = 1;
1103         } else {
1104             alg = 2;
1105         }
1106     } else if (communicator_size < 32) {
1107         if (per_rank_dsize < 128) {
1108             alg = 1;
1109         } else if (per_rank_dsize < 262144) {
1110             alg = 2;
1111         } else {
1112             alg = 3;
1113         }
1114     } else if (communicator_size < 64) {
1115         if (per_rank_dsize < 256) {
1116             alg = 1;
1117         } else if (per_rank_dsize < 8192) {
1118             alg = 2;
1119         } else {
1120             alg = 3;
1121         }
1122     } else if (communicator_size < 128) {
1123         if (per_rank_dsize < 256) {
1124             alg = 1;
1125         } else if (per_rank_dsize < 4096) {
1126             alg = 2;
1127         } else {
1128             alg = 3;
1129         }
1130     } else if (communicator_size < 256) {
1131         if (per_rank_dsize < 1024) {
1132             alg = 2;
1133         } else if (per_rank_dsize < 65536) {
1134             alg = 4;
1135         } else {
1136             alg = 3;
1137         }
1138     } else if (communicator_size < 512) {
1139         if (per_rank_dsize < 1024) {
1140             alg = 2;
1141         } else {
1142             alg = 3;
1143         }
1144     } else if (communicator_size < 1024) {
1145         if (per_rank_dsize < 512) {
1146             alg = 2;
1147         } else if (per_rank_dsize < 1024) {
1148             alg = 1;
1149         } else if (per_rank_dsize < 4096) {
1150             alg = 2;
1151         } else if (per_rank_dsize < 1048576) {
1152             alg = 4;
1153         } else {
1154             alg = 3;
1155         }
1156     } else {
1157         if (per_rank_dsize < 4096) {
1158             alg = 2;
1159         } else {
1160             alg = 4;
1161         }
1162     }
1163
1164     return funcs[alg-1](sbuf, scount, sdtype, rbuf, rcounts, rdispls, rdtype, comm);
1165 }
1166
1167 int gather__ompi(const void *sbuf, int scount,
1168                  MPI_Datatype sdtype,
1169                  void* rbuf, int rcount,
1170                  MPI_Datatype rdtype,
1171                  int root,
1172                  MPI_Comm  comm
1173                  )
1174 {
1175     int communicator_size, rank;
1176     size_t dsize, total_dsize;
1177     int alg = 1;
1178     communicator_size = comm->size();
1179     rank = comm->rank();
1180
1181     if (rank == root) {
1182         dsize = rdtype->size();
1183         total_dsize = dsize * rcount;
1184     } else {
1185         dsize = sdtype->size();
1186         total_dsize = dsize * scount;
1187     }
1188     int (*funcs[])(const void*, int, MPI_Datatype, void*, int, MPI_Datatype, int, MPI_Comm) = {
1189         &gather__ompi_basic_linear,
1190         &gather__ompi_binomial,
1191         &gather__ompi_linear_sync
1192     };
1193     /** Algorithms:
1194      *  {1, "basic_linear"},
1195      *  {2, "binomial"},
1196      *  {3, "linear_sync"},
1197      *
1198      * We do not make any rank specific checks since the params
1199      * should be uniform across ranks.
1200      */
1201     if (communicator_size < 4) {
1202         if (total_dsize < 2) {
1203             alg = 3;
1204         } else if (total_dsize < 4) {
1205             alg = 1;
1206         } else if (total_dsize < 32768) {
1207             alg = 2;
1208         } else if (total_dsize < 65536) {
1209             alg = 1;
1210         } else if (total_dsize < 131072) {
1211             alg = 2;
1212         } else {
1213             alg = 3;
1214         }
1215     } else if (communicator_size < 8) {
1216         if (total_dsize < 1024) {
1217             alg = 2;
1218         } else if (total_dsize < 8192) {
1219             alg = 1;
1220         } else if (total_dsize < 32768) {
1221             alg = 2;
1222         } else if (total_dsize < 262144) {
1223             alg = 1;
1224         } else {
1225             alg = 3;
1226         }
1227     } else if (communicator_size < 256) {
1228         alg = 2;
1229     } else if (communicator_size < 512) {
1230         if (total_dsize < 2048) {
1231             alg = 2;
1232         } else if (total_dsize < 8192) {
1233             alg = 1;
1234         } else {
1235             alg = 2;
1236         }
1237     } else {
1238         alg = 2;
1239     }
1240
1241     return funcs[alg-1](sbuf, scount, sdtype, rbuf, rcount, rdtype, root, comm);
1242 }
1243
1244
1245 int scatter__ompi(const void *sbuf, int scount,
1246                   MPI_Datatype sdtype,
1247                   void* rbuf, int rcount,
1248                   MPI_Datatype rdtype,
1249                   int root, MPI_Comm  comm
1250                   )
1251 {
1252     int communicator_size, rank;
1253     size_t dsize, total_dsize;
1254     int alg = 1;
1255
1256     communicator_size = comm->size();
1257     rank = comm->rank();
1258     if (root == rank) {
1259         dsize=sdtype->size();
1260         total_dsize = dsize * scount;
1261     } else {
1262         dsize=rdtype->size();
1263         total_dsize = dsize * rcount;
1264     }
1265     int (*funcs[])(const void*, int, MPI_Datatype, void*, int, MPI_Datatype, int, MPI_Comm) = {
1266         &scatter__ompi_basic_linear,
1267         &scatter__ompi_binomial,
1268         &scatter__ompi_linear_nb
1269     };
1270     /** Algorithms:
1271      *  {1, "basic_linear"},
1272      *  {2, "binomial"},
1273      *  {3, "linear_nb"},
1274      *
1275      * We do not make any rank specific checks since the params
1276      * should be uniform across ranks.
1277      */
1278     if (communicator_size < 4) {
1279         if (total_dsize < 2) {
1280             alg = 3;
1281         } else if (total_dsize < 131072) {
1282             alg = 1;
1283         } else if (total_dsize < 262144) {
1284             alg = 3;
1285         } else {
1286             alg = 1;
1287         }
1288     } else if (communicator_size < 8) {
1289         if (total_dsize < 2048) {
1290             alg = 2;
1291         } else if (total_dsize < 4096) {
1292             alg = 1;
1293         } else if (total_dsize < 8192) {
1294             alg = 2;
1295         } else if (total_dsize < 32768) {
1296             alg = 1;
1297         } else if (total_dsize < 1048576) {
1298             alg = 3;
1299         } else {
1300             alg = 1;
1301         }
1302     } else if (communicator_size < 16) {
1303         if (total_dsize < 16384) {
1304             alg = 2;
1305         } else if (total_dsize < 1048576) {
1306             alg = 3;
1307         } else {
1308             alg = 1;
1309         }
1310     } else if (communicator_size < 32) {
1311         if (total_dsize < 16384) {
1312             alg = 2;
1313         } else if (total_dsize < 32768) {
1314             alg = 1;
1315         } else {
1316             alg = 3;
1317         }
1318     } else if (communicator_size < 64) {
1319         if (total_dsize < 512) {
1320             alg = 2;
1321         } else if (total_dsize < 8192) {
1322             alg = 3;
1323         } else if (total_dsize < 16384) {
1324             alg = 2;
1325         } else {
1326             alg = 3;
1327         }
1328     } else {
1329         if (total_dsize < 512) {
1330             alg = 2;
1331         } else {
1332             alg = 3;
1333         }
1334     }
1335
1336     return funcs[alg-1](sbuf, scount, sdtype, rbuf, rcount, rdtype, root, comm);
1337 }
1338
1339 }
1340 }