Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
6bec031ed8b206dbc0a7e0317c2a2dfcdca1cad6
[simgrid.git] / src / smpi / colls / smpi_intel_mpi_selector.cpp
1 /* selector for collective algorithms based on openmpi's default coll_tuned_decision_fixed selector */
2
3 /* Copyright (c) 2009-2010, 2013-2017. The SimGrid Team.
4  * All rights reserved.                                                     */
5
6 /* This program is free software; you can redistribute it and/or modify it
7  * under the terms of the license (GNU LGPL) which comes with this package. */
8
9 #include "colls_private.h"
10
11 // This selector is based on information gathered on the Stampede cluster, with Intel MPI 4.1.3.049, and from the intel reference manual. The data was gathered launching runs with 1,2,4,8,16 processes per node.
12
13 #define INTEL_MAX_NB_THRESHOLDS  32
14 #define INTEL_MAX_NB_NUMPROCS  12
15 #define INTEL_MAX_NB_PPN  5  /* 1 2 4 8 16 ppn */
16
17 typedef struct {
18   unsigned int max_size;
19   int algo;
20 } intel_tuning_table_size_element;
21
22 typedef struct {
23   int max_num_proc;
24   int num_elems;
25   intel_tuning_table_size_element elems[INTEL_MAX_NB_THRESHOLDS];
26 } intel_tuning_table_numproc_element;
27
28 typedef struct {
29   int ppn;
30   intel_tuning_table_numproc_element elems[INTEL_MAX_NB_NUMPROCS];
31 } intel_tuning_table_element;
32
33 /*
34 I_MPI_ADJUST_ALLREDUCE
35
36 MPI_Allreduce
37
38 1 - Recursive doubling algorithm
39 2 - Rabenseifner's algorithm
40 3 - Reduce + Bcast algorithm
41 4 - Topology aware Reduce + Bcast algorithm
42 5 - Binomial gather + scatter algorithm
43 6 - Topology aware binominal gather + scatter algorithm
44 7 - Shumilin's ring algorithm
45 8 - Ring algorithm
46
47   as Shumilin's ring algorithm is unknown, default to ring'
48 */
49
50 namespace simgrid{
51 namespace smpi{
52
53 int (*intel_allreduce_functions_table[])(void *sendbuf,
54       void *recvbuf,
55       int count,
56       MPI_Datatype datatype,
57       MPI_Op op, MPI_Comm comm) ={
58       Coll_allreduce_rdb::allreduce,
59       Coll_allreduce_rab1::allreduce,
60       Coll_allreduce_redbcast::allreduce,
61       Coll_allreduce_mvapich2_two_level::allreduce,
62       Coll_allreduce_smp_binomial::allreduce,
63       Coll_allreduce_mvapich2_two_level::allreduce,
64       Coll_allreduce_ompi_ring_segmented::allreduce,
65       Coll_allreduce_ompi_ring_segmented::allreduce
66 };
67
68 intel_tuning_table_element intel_allreduce_table[] =
69 {
70   {1,{
71     { 2,9,{
72       {6,7},
73       {85,1},
74       {192,7},
75       {853,1},
76       {1279,7},
77       {16684,1},
78       {34279,8},
79       {1681224,3},
80       {2147483647,7}
81     }
82     },
83     { 4, 8,{
84       {16,7},
85       {47,1},
86       {2062,7},
87       {16699,1},
88       {33627,7},
89       {70732,8},
90       {1300705,3},
91       {2147483647,8}
92     }
93     },
94     {8,8,{
95       {118,1},
96       {146,4},
97       {16760,1},
98       {36364,6},
99       {136239,8},
100       {315710,7},
101       {3220366,3},
102       {2147483647,8}
103       }
104     },
105     {16,7,{
106       {934,1},
107       {1160,6},
108       {15505,1},
109       {52730,2},
110       {300705,8},
111       {563680,7},
112       {2147483647,3}
113       }
114     },
115     {2147483647,11,{
116       {5,6},
117       {11,4},
118       {182,1},
119       {700,6},
120       {1450,4},
121       {11146,1},
122       {25539,6},
123       {37634,4},
124       {93784,6},
125       {817658,2},
126       {2147483647,3}
127     }
128     }
129   }
130   },
131   {2,{
132     { 4,6,{
133       {2084,7},
134       {15216,1},
135       {99715,7},
136       {168666,3},
137       {363889,2},
138       {2147483647,7}
139     }
140     },
141     { 8,6,{
142       {14978,1},
143       {66879,2},
144       {179296,8},
145       {304801,3},
146       {704509,7},
147       {2147483647,2}
148     }
149     },
150     { 16,6,{
151       {16405,1},
152       {81784,2},
153       {346385,8},
154       {807546,7},
155       {1259854,2},
156       {2147483647,3}
157     }
158     },
159     { 32,4,{
160       {8913,1},
161       {103578,2},
162       {615876,8},
163       {2147483647,2}
164     }
165     },
166     { 64,7,{
167       {1000,1},
168       {2249,2},
169       {6029,1},
170       {325357,2},
171       {1470976,8},
172       {2556670,7},
173       {2147483647,3}
174     }
175     },
176     { 128,5,{
177       {664,1},
178       {754706,2},
179       {1663862,4},
180       {3269097,2},
181       {2147483647,7}
182     }
183     },
184     { 2147483647,3,{
185       {789,1},
186       {2247589,2},
187       {2147483647,8}
188     }
189     }
190   }
191   },
192   {4,{
193     { 4,4,{
194       {5738,1},
195       {197433,2},
196       {593742,7},
197       {2147483647,2}
198     }
199     },
200     { 8,7,{
201       {5655,1},
202       {75166,2},
203       {177639,8},
204       {988014,3},
205       {1643869,2},
206       {2494859,8},
207       {2147483647,2}
208     }
209     },
210     { 16,7,{
211       {587,1},
212       {3941,2},
213       {9003,1},
214       {101469,2},
215       {355768,8},
216       {3341814,3},
217       {2147483647,8}
218     }
219     },
220     { 32,4,{
221       {795,1},
222       {146567,2},
223       {732118,8},
224       {2147483647,3}
225     }
226     },
227     { 64,4,{
228       {528,1},
229       {221277,2},
230       {1440737,8},
231       {2147483647,3}
232     }
233     },
234     { 128,4,{
235       {481,1},
236       {593833,2},
237       {2962021,8},
238       {2147483647,7}
239     }
240     },
241     { 256,2,{
242       {584,1},
243       {2147483647,2}
244     }
245     },
246     { 2147483647,3,{
247       {604,1},
248       {2997006,2},
249       {2147483647,8}
250     }
251     }
252   }
253   },
254   {8,{
255     { 8,6,{
256       {2560,1},
257       {114230,6},
258       {288510,8},
259       {664038,2},
260       {1339913,6},
261       {2147483647,4}
262     }
263     },
264     { 16,5,{
265       {497,1},
266       {54201,2},
267       {356217,8},
268       {3413609,3},
269       {2147483647,8}
270     }
271     },
272     { 32,5,{
273       {377,1},
274       {109745,2},
275       {716514,8},
276       {3976768,3},
277       {2147483647,8}
278     }
279     },
280     { 64,6,{
281       {109,1},
282       {649,5},
283       {266080,2},
284       {1493331,8},
285       {2541403,7},
286       {2147483647,3}
287     }
288     },
289     { 128,4,{
290       {7,1},
291       {751,5},
292       {408808,2},
293       {2147483647,8}
294     }
295     },
296     { 256,3,{
297       {828,5},
298       {909676,2},
299       {2147483647,8}
300     }
301     },
302     { 512,5,{
303       {847,5},
304       {1007066,2},
305       {1068775,4},
306       {2803389,2},
307       {2147483647,8}
308     }
309     },
310     { 2147483647,3,{
311       {1974,5},
312       {4007876,2},
313       {2147483647,8}
314     }
315     }
316   }
317   },
318   {16,{
319     { 16,12,{
320       {409,1},
321       {768,6},
322       {1365,4},
323       {3071,6},
324       {11299,2},
325       {21746,6},
326       {55629,2},
327       {86065,4},
328       {153867,2},
329       {590560,6},
330       {1448760,2},
331       {2147483647,8},
332     }
333     },
334     { 32,8,{
335       {6,1},
336       {24,5},
337       {86,1},
338       {875,5},
339       {74528,2},
340       {813050,8},
341       {1725981,7},
342       {2147483647,8},
343     }
344     },
345     { 64,6,{
346       {1018,5},
347       {1217,6},
348       {2370,5},
349       {160654,2},
350       {1885487,8},
351       {2147483647,3},
352     }
353     },
354     { 128,4,{
355       {2291,5},
356       {434465,2},
357       {3525103,8},
358       {2147483647,7},
359     }
360     },
361     { 256,3,{
362       {2189,5},
363       {713154,2},
364       {2147483647,8},
365     }
366     },
367     { 512,3,{
368       {2140,5},
369       {1235056,2},
370       {2147483647,8},
371     }
372     },
373     { 2147483647,3,{
374       {2153,5},
375       {2629855,2},
376       {2147483647,8},
377     }
378     }
379   }
380   }
381 };
382
383
384
385 /*I_MPI_ADJUST_ALLTOALL
386
387 MPI_Alltoall
388
389 1. Bruck's algorithm
390 2. Isend/Irecv + waitall algorithm
391 3. Pair wise exchange algorithm
392 4. Plum's algorithm
393
394 */
395
396
397 intel_tuning_table_element intel_alltoall_table[] =
398 {
399   {1,{
400     { 2,1,
401         {
402         {2147483647,3}
403         }
404     },
405     { 4,2,
406         {
407         {0,4},
408         {2147483647,2}
409         }
410     },
411     {8,1,
412         {
413         {2147483647,2}
414         }
415     },
416     {16,5,
417         {
418         {0,3},
419         {84645,2},
420         {167570,3},
421         {413152,4},
422         {2147483647,2}
423         }
424     },
425     {32,6,
426         {
427         {61,1},
428         {164,2},
429         {696,1},
430         {143254,2},
431         {387024,3},
432         {2147483647,2}
433         },
434     },
435     {64,4,
436         {
437         {523,1},
438         {146088,2},
439         {488989,4},
440         {2147483647,2}
441         }
442     },
443     {2147483647,3,
444         {
445         {270,1},
446         {628,4},
447         {2147483647,2}
448         }
449     }
450   }
451   },
452   {2, {
453     { 4,4,{
454       {1,2},
455       {75,3},
456       {131072,2},
457       {2147483647,2}
458     }
459     },
460     { 8,3,{
461       {709,1},
462       {131072,2},
463       {2147483647,2}
464     }
465     },
466     { 16,4,{
467       {40048,2},
468       {131072,3},
469       {155927,3},
470       {2147483647,4}
471     }
472     },
473     { 32,7,{
474       {105,1},
475       {130,2},
476       {1030,1},
477       {58893,2},
478       {131072,2},
479       {271838,3},
480       {2147483647,2}
481     }
482     },
483     { 2147483647,8,{
484       {521,1},
485       {2032,4},
486       {2412,2},
487       {4112,4},
488       {61620,2},
489       {131072,3},
490       {427408,3},
491       {2147483647,4}
492     }
493     }
494   }
495   },
496   {4,{
497     { 8,3,{
498       {512,1},
499       {32768,2},
500       {2147483647,2}
501     }
502     },
503     { 64,8,{
504       {7,1},
505       {199,4},
506       {764,1},
507       {6409,4},
508       {20026,2},
509       {32768,3},
510       {221643,4},
511       {2147483647,3}
512     }
513     },
514     { 2147483647,7,{
515       {262,1},
516       {7592,4},
517       {22871,2},
518       {32768,3},
519       {47538,3},
520       {101559,4},
521       {2147483647,3}
522     }
523     }
524   }
525   },
526   {8,{
527     { 16,6,{
528       {973,1},
529       {5126,4},
530       {16898,2},
531       {32768,4},
532       {65456,4},
533       {2147483647,2}
534     }
535     },
536     { 32,7,{
537       {874,1},
538       {6727,4},
539       {17912,2},
540       {32768,3},
541       {41513,3},
542       {199604,4},
543       {2147483647,3}
544     }
545     },
546     { 64,8,{
547       {5,1},
548       {114,4},
549       {552,1},
550       {8130,4},
551       {32768,3},
552       {34486,3},
553       {160113,4},
554       {2147483647,3}
555     }
556     },
557     { 128,6,{
558       {270,1},
559       {3679,4},
560       {32768,3},
561       {64367,3},
562       {146595,4},
563       {2147483647,3}
564     }
565     },
566     { 2147483647,4,{
567       {133,1},
568       {4017,4},
569       {32768,3},
570       {76351,4},
571       {2147483647,3}
572     }
573     }
574   }
575   },
576   {16,{
577     { 32,7,{
578       {963,1},
579       {1818,4},
580       {20007,2},
581       {32768,4},
582       {54296,4},
583       {169735,3},
584       {2147483647,2}
585     }
586     },
587     { 64,11,{
588       {17,1},
589       {42,4},
590       {592,1},
591       {2015,4},
592       {2753,2},
593       {6496,3},
594       {20402,4},
595       {32768,3},
596       {36246,3},
597       {93229,4},
598       {2147483647,3}
599     }
600     },
601     { 128,9,{
602       {18,1},
603       {40,4},
604       {287,1},
605       {1308,4},
606       {6842,1},
607       {32768,3},
608       {36986,3},
609       {129081,4},
610       {2147483647,3}
611     }
612     },
613     { 256,7,{
614       {135,1},
615       {1538,4},
616       {3267,1},
617       {4132,3},
618       {31469,4},
619       {32768,3},
620       {2147483647,3}
621     }
622     },
623     { 2147483647,8,{
624       {66,1},
625       {1637,4},
626       {2626,1},
627       {4842,4},
628       {32768,3},
629       {33963,3},
630       {72978,4},
631       {2147483647,3}
632     }
633     }
634   }
635   }
636 };
637 int (*intel_alltoall_functions_table[])(void *sbuf, int scount,
638                                              MPI_Datatype sdtype,
639                                              void* rbuf, int rcount,
640                                              MPI_Datatype rdtype,
641                                              MPI_Comm comm) ={
642       Coll_alltoall_bruck::alltoall,
643       Coll_alltoall_mvapich2_scatter_dest::alltoall,
644       Coll_alltoall_pair::alltoall,
645       Coll_alltoall_mvapich2::alltoall//Plum is proprietary ? (and super efficient)
646 };
647
648 /*I_MPI_ADJUST_BARRIER
649
650 MPI_Barrier
651
652 1. Dissemination algorithm
653 2. Recursive doubling algorithm
654 3. Topology aware dissemination algorithm
655 4. Topology aware recursive doubling algorithm
656 5. Binominal gather + scatter algorithm
657 6. Topology aware binominal gather + scatter algorithm
658
659 */
660 static int intel_barrier_gather_scatter(MPI_Comm comm){
661     //our default barrier performs a antibcast/bcast
662     Coll_barrier_default::barrier(comm);
663     return MPI_SUCCESS;
664 }
665
666 int (*intel_barrier_functions_table[])(MPI_Comm comm) ={
667       Coll_barrier_ompi_basic_linear::barrier,
668       Coll_barrier_ompi_recursivedoubling::barrier,
669       Coll_barrier_ompi_basic_linear::barrier,
670       Coll_barrier_ompi_recursivedoubling::barrier,
671       intel_barrier_gather_scatter,
672       intel_barrier_gather_scatter
673 };
674
675 intel_tuning_table_element intel_barrier_table[] =
676 {
677   {1,{
678     {2,1,
679         {
680         {2147483647,2}
681         }
682     },
683     {4,1,
684         {
685         {2147483647,6}
686         }
687     },
688     {8,1,
689         {
690         {2147483647,1}
691         }
692     },
693     {64,1,
694         {
695         {2147483647,2}
696         }
697     },
698     {2147483647,1,
699         {
700         {2147483647,6}
701         }
702     }
703   }
704   },
705   {2,{
706     { 2,1,{
707       {2147483647,1}
708     }
709     },
710     { 4,1,{
711       {2147483647,3}
712     }
713     },
714     { 8,1,{
715       {2147483647,5}
716     }
717     },
718     { 32,1,{
719       {2147483647,2}
720     }
721     },
722     { 128,1,{
723       {2147483647,3}
724     }
725     },
726     { 2147483647,1,{
727       {2147483647,4}
728     }
729     }
730   }
731   },
732   {4,{
733     { 4,1,{
734       {2147483647,2}
735     }
736     },
737     { 8,1,{
738       {2147483647,5}
739     }
740     },
741     { 32,1,{
742       {2147483647,2}
743     }
744     },
745     { 2147483647,1,{
746       {2147483647,4}
747     }
748     }
749   }
750   },
751   {8,{
752     { 8,1,{
753       {2147483647,1}
754     }
755     },
756     { 32,1,{
757       {2147483647,2}
758     }
759     },
760     { 2147483647,1,{
761       {2147483647,4}
762     }
763     }
764   }
765   },
766   {16,{
767     { 4,1,{
768       {2147483647,2}
769     }
770     },
771     { 8,1,{
772       {2147483647,5}
773     }
774     },
775     { 32,1,{
776       {2147483647,2}
777     }
778     },
779     { 2147483647,1,{
780       {2147483647,4}
781     }
782     }
783   }
784   }
785 };
786
787
788 /*I_MPI_ADJUST_BCAST
789
790 MPI_Bcast
791
792 1. Binomial algorithm
793 2. Recursive doubling algorithm
794 3. Ring algorithm
795 4. Topology aware binomial algorithm
796 5. Topology aware recursive doubling algorithm
797 6. Topology aware ring algorithm
798 7. Shumilin's bcast algorithm
799 */
800
801 int (*intel_bcast_functions_table[])(void *buff, int count,
802                                           MPI_Datatype datatype, int root,
803                                           MPI_Comm  comm) ={
804       Coll_bcast_binomial_tree::bcast,
805       //Coll_bcast_scatter_rdb_allgather::bcast,
806       Coll_bcast_NTSL::bcast,
807       Coll_bcast_NTSL::bcast,
808       Coll_bcast_SMP_binomial::bcast,
809       //Coll_bcast_scatter_rdb_allgather::bcast,
810       Coll_bcast_NTSL::bcast,
811       Coll_bcast_SMP_linear::bcast,
812       Coll_bcast_mvapich2::bcast,//we don't know shumilin's algo'
813 };
814
815 intel_tuning_table_element intel_bcast_table[] =
816 {
817   {1,{
818     {2,9,
819         {
820         {1,2},
821         {402,7},
822         {682,5},
823         {1433,4},
824         {5734,7},
825         {21845,1},
826         {95963,6},
827         {409897,5},
828         {2147483647,1}
829         }
830     },
831     {4,1,
832         {
833         {2147483647,7}
834         }
835     },
836     {8,11,
837         {
838         {3,6},
839         {4,7},
840         {25,6},
841         {256,1},
842         {682,6},
843         {1264,1},
844         {2234,6},
845         {6655,5},
846         {16336,1},
847         {3998434,7},
848         {2147483647,6}
849         }
850     },
851     {2147483647,1,
852         {
853         {2147483647,7}
854         }
855     }
856   }
857   },
858   {2,{
859     { 4,6,{
860       {806,4},
861       {18093,7},
862       {51366,6},
863       {182526,4},
864       {618390,1},
865       {2147483647,7}
866     }
867     },
868     { 8,6,{
869       {24,1},
870       {74,4},
871       {18137,1},
872       {614661,7},
873       {1284626,1},
874       {2147483647,2}
875     }
876     },
877     { 16,4,{
878       {1,1},
879       {158,7},
880       {16955,1},
881       {2147483647,7}
882     }
883     },
884     { 32,3,{
885       {242,7},
886       {10345,1},
887       {2147483647,7}
888     }
889     },
890     { 2147483647,4,{
891       {1,1},
892       {737,7},
893       {5340,1},
894       {2147483647,7}
895     }
896     }
897   }
898   },
899   {4,{
900     { 8,4,{
901       {256,4},
902       {17181,1},
903       {1048576,7},
904       {2147483647,7}
905     }
906     },
907     { 2147483647,1,{
908       {2147483647,7}
909     }
910     }
911   }
912   },
913   {8,{
914     { 16,5,{
915       {3,1},
916       {318,7},
917       {1505,1},
918       {1048576,7},
919       {2147483647,7}
920     }
921     },
922     { 32,3,{
923       {422,7},
924       {851,1},
925       {2147483647,7}
926     }
927     },
928     { 64,3,{
929       {468,7},
930       {699,1},
931       {2147483647,7}
932     }
933     },
934     { 2147483647,1,{
935       {2147483647,7}
936     }
937     }
938   }
939   },
940   {16,{
941     { 8,4,{
942       {256,4},
943       {17181,1},
944       {1048576,7},
945       {2147483647,7}
946     }
947     },
948     { 2147483647,1,{
949       {2147483647,7}
950     }
951     }
952   }
953   }
954 };
955
956
957 /*I_MPI_ADJUST_REDUCE
958
959 MPI_Reduce
960
961 1. Shumilin's algorithm
962 2. Binomial algorithm
963 3. Topology aware Shumilin's algorithm
964 4. Topology aware binomial algorithm
965 5. Rabenseifner's algorithm
966 6. Topology aware Rabenseifner's algorithm
967
968 */
969
970 int (*intel_reduce_functions_table[])(void *sendbuf, void *recvbuf,
971                                             int count, MPI_Datatype  datatype,
972                                             MPI_Op   op, int root,
973                                             MPI_Comm   comm) ={
974       Coll_reduce_mvapich2::reduce,
975       Coll_reduce_binomial::reduce,
976       Coll_reduce_mvapich2::reduce,
977       Coll_reduce_mvapich2_two_level::reduce,
978       Coll_reduce_rab::reduce,
979       Coll_reduce_rab::reduce
980 };
981
982 intel_tuning_table_element intel_reduce_table[] =
983 {
984   {1,{
985     {2147483647,1,
986       {
987       {2147483647,1}
988       }
989     }
990   }
991   },
992   {2,{
993     { 2,1,{
994       {2147483647,1}
995     }
996     },
997     { 4,2,{
998       {10541,3},
999       {2147483647,1}
1000     }
1001     },
1002     { 2147483647,1,{
1003       {2147483647,1}
1004     }
1005     }
1006   }
1007   },
1008   {4,{
1009     { 256,1,{
1010       {2147483647,1}
1011     }
1012     },
1013     { 2147483647,2,{
1014       {45,3},
1015       {2147483647,1}
1016     }
1017     }
1018   }
1019   },
1020   {8,{
1021     { 512,1,{
1022       {2147483647,1}
1023     }
1024     },
1025     { 2147483647,3,{
1026       {5,1},
1027       {11882,3},
1028       {2147483647,1}
1029     }
1030     }
1031   }
1032   },
1033   {16,{
1034     { 256,1,{
1035       {2147483647,1}
1036     }
1037     },
1038     { 2147483647,2,{
1039       {45,3},
1040       {2147483647,1}
1041     }
1042     }
1043   }
1044   }
1045 };
1046
1047 /* I_MPI_ADJUST_REDUCE_SCATTER
1048
1049 MPI_Reduce_scatter
1050
1051 1. Recursive having algorithm
1052 2. Pair wise exchange algorithm
1053 3. Recursive doubling algorithm
1054 4. Reduce + Scatterv algorithm
1055 5. Topology aware Reduce + Scatterv algorithm
1056
1057 */
1058 static  int intel_reduce_scatter_reduce_scatterv(void *sbuf, void *rbuf,
1059                                                     int *rcounts,
1060                                                     MPI_Datatype dtype,
1061                                                     MPI_Op  op,
1062                                                     MPI_Comm  comm)
1063 {
1064   Coll_reduce_scatter_default::reduce_scatter(sbuf, rbuf, rcounts,dtype, op,comm);
1065   return MPI_SUCCESS;
1066 }
1067
1068 static  int  intel_reduce_scatter_recursivehalving(void *sbuf, void *rbuf,
1069                                                     int *rcounts,
1070                                                     MPI_Datatype dtype,
1071                                                     MPI_Op  op,
1072                                                     MPI_Comm  comm)
1073 {
1074   if(op==MPI_OP_NULL || op->is_commutative())
1075     return Coll_reduce_scatter_ompi_basic_recursivehalving::reduce_scatter(sbuf, rbuf, rcounts,dtype, op,comm);
1076   else
1077     return Coll_reduce_scatter_mvapich2::reduce_scatter(sbuf, rbuf, rcounts,dtype, op,comm);
1078 }
1079
1080 int (*intel_reduce_scatter_functions_table[])( void *sbuf, void *rbuf,
1081                                                     int *rcounts,
1082                                                     MPI_Datatype dtype,
1083                                                     MPI_Op  op,
1084                                                     MPI_Comm  comm
1085                                                     ) ={
1086       intel_reduce_scatter_recursivehalving,
1087       Coll_reduce_scatter_mpich_pair::reduce_scatter,
1088       Coll_reduce_scatter_mpich_rdb::reduce_scatter,
1089       intel_reduce_scatter_reduce_scatterv,
1090       intel_reduce_scatter_reduce_scatterv
1091 };
1092
1093 intel_tuning_table_element intel_reduce_scatter_table[] =
1094 {
1095   {1,{
1096     {2,5,
1097     {
1098         {5,4},
1099         {522429,2},
1100         {1375877,5},
1101         {2932736,2},
1102         {2147483647,5}
1103         }
1104     },
1105     {4,9,
1106         {
1107         {4,4},
1108         {15,1},
1109         {120,3},
1110         {651,1},
1111         {12188,3},
1112         {33890,1},
1113         {572117,2},
1114         {1410202,5},
1115         {2147483647,2}
1116         }
1117     },
1118     {8,7,
1119         {
1120         {4,4},
1121         {2263,1},
1122         {25007,3},
1123         {34861,1},
1124         {169625,2},
1125         {2734000,4},
1126         {2147483647,2}
1127         }
1128     },
1129     {16,5,
1130         {
1131         {4,4},
1132         {14228,1},
1133         {46084,3},
1134         {522139,2},
1135         {2147483647,5}
1136         }
1137     },
1138     {32,5,
1139         {
1140         {4,4},
1141         {27516,1},
1142         {61693,3},
1143         {2483469,2},
1144         {2147483647,5}
1145         }
1146     },
1147     {64,4,
1148         {
1149         {0,3},
1150         {4,4},
1151         {100396,1},
1152         {2147483647,2}
1153         }
1154     },
1155     {2147483647,6,
1156         {
1157         {0,3},
1158         {4,4},
1159         {186926,1},
1160         {278259,3},
1161         {1500100,2},
1162         {2147483647,5}
1163         }
1164     }
1165   }
1166   },
1167   {2,{
1168     { 2,2,{
1169       {6,1},
1170       {2147483647,2}
1171     }
1172     },
1173     { 4,7,{
1174       {5,4},
1175       {13,5},
1176       {59,3},
1177       {76,1},
1178       {91488,3},
1179       {680063,4},
1180       {2147483647,2}
1181     }
1182     },
1183     { 8,8,{
1184       {4,4},
1185       {11,5},
1186       {31,1},
1187       {69615,3},
1188       {202632,2},
1189       {396082,5},
1190       {1495696,4},
1191       {2147483647,2}
1192     }
1193     },
1194     { 16,1,{
1195       {4,4},
1196       {345,1},
1197       {79523,3},
1198       {2147483647,2}
1199     }
1200     },
1201     { 32,5,{
1202       {0,3},
1203       {4,4},
1204       {992,1},
1205       {71417,3},
1206       {2147483647,2}
1207     }
1208     },
1209     { 64,4,{
1210       {4,4},
1211       {1472,1},
1212       {196592,3},
1213       {2147483647,2}
1214     }
1215     },
1216     { 128,5,{
1217       {0,3},
1218       {4,4},
1219       {32892,1},
1220       {381072,3},
1221       {2147483647,2}
1222     }
1223     },
1224     { 2147483647,6,{
1225       {0,2},
1226       {4,4},
1227       {33262,1},
1228       {1571397,3},
1229       {2211398,5},
1230       {2147483647,4}
1231     }
1232     }
1233   }
1234   },
1235   {4,{
1236     { 4,7,{
1237       {12,4},
1238       {27,5},
1239       {49,3},
1240       {187,1},
1241       {405673,3},
1242       {594687,4},
1243       {2147483647,2}
1244     }
1245     },
1246     { 8,5,{
1247       {24,5},
1248       {155,1},
1249       {204501,3},
1250       {274267,5},
1251       {2147483647,4}
1252     }
1253     },
1254     { 16,6,{
1255       {63,1},
1256       {72,3},
1257       {264,1},
1258       {168421,3},
1259       {168421,4},
1260       {2147483647,2}
1261     }
1262     },
1263     { 32,10,{
1264       {0,3},
1265       {4,4},
1266       {12,1},
1267       {18,5},
1268       {419,1},
1269       {188739,3},
1270       {716329,4},
1271       {1365841,5},
1272       {2430194,2},
1273       {2147483647,4}
1274     }
1275     },
1276     { 64,8,{
1277       {0,3},
1278       {4,4},
1279       {17,5},
1280       {635,1},
1281       {202937,3},
1282       {308253,5},
1283       {1389874,4},
1284       {2147483647,2}
1285     }
1286     },
1287     { 128,8,{
1288       {0,3},
1289       {4,4},
1290       {16,5},
1291       {1238,1},
1292       {280097,3},
1293       {631434,5},
1294       {2605072,4},
1295       {2147483647,2}
1296     }
1297     },
1298     { 256,5,{
1299       {0,2},
1300       {4,4},
1301       {16,5},
1302       {2418,1},
1303       {2147483647,3}
1304     }
1305     },
1306     { 2147483647,6,{
1307       {0,2},
1308       {4,4},
1309       {16,5},
1310       {33182,1},
1311       {3763779,3},
1312       {2147483647,4}
1313     }
1314     }
1315   }
1316   },
1317   {8,{
1318     { 8,6,{
1319       {5,4},
1320       {494,1},
1321       {97739,3},
1322       {522836,2},
1323       {554174,5},
1324       {2147483647,2}
1325     }
1326     },
1327     { 16,8,{
1328       {5,4},
1329       {62,1},
1330       {94,3},
1331       {215,1},
1332       {185095,3},
1333       {454784,4},
1334       {607911,5},
1335       {2147483647,4}
1336     }
1337     },
1338     { 32,7,{
1339       {0,3},
1340       {4,4},
1341       {302,1},
1342       {250841,3},
1343       {665822,4},
1344       {1760980,5},
1345       {2147483647,4}
1346     }
1347     },
1348     { 64,8,{
1349       {0,3},
1350       {4,4},
1351       {41,5},
1352       {306,1},
1353       {332405,3},
1354       {1269189,4},
1355       {3712421,5},
1356       {2147483647,4}
1357     }
1358     },
1359     { 128,6,{
1360       {0,3},
1361       {4,4},
1362       {39,5},
1363       {526,1},
1364       {487878,3},
1365       {2147483647,4}
1366     }
1367     },
1368     { 256,8,{
1369       {0,2},
1370       {4,4},
1371       {36,5},
1372       {1382,1},
1373       {424162,3},
1374       {632881,5},
1375       {1127566,3},
1376       {2147483647,4}
1377     }
1378     },
1379     { 512,4,{
1380       {4,4},
1381       {34,5},
1382       {5884,1},
1383       {2147483647,3}
1384     }
1385     },
1386     { 2147483647,4,{
1387       {5,4},
1388       {32,5},
1389       {25105,1},
1390       {2147483647,3}
1391     }
1392     }
1393   }
1394   },
1395   {16,{
1396     { 4,7,{
1397       {12,4},
1398       {27,5},
1399       {49,3},
1400       {187,1},
1401       {405673,3},
1402       {594687,4},
1403       {2147483647,2}
1404     }
1405     },
1406     { 8,5,{
1407       {24,5},
1408       {155,1},
1409       {204501,3},
1410       {274267,5},
1411       {2147483647,4}
1412     }
1413     },
1414     { 16,6,{
1415       {63,1},
1416       {72,3},
1417       {264,1},
1418       {168421,3},
1419       {168421,4},
1420       {2147483647,2}
1421     }
1422     },
1423     { 32,10,{
1424       {0,3},
1425       {4,4},
1426       {12,1},
1427       {18,5},
1428       {419,1},
1429       {188739,3},
1430       {716329,4},
1431       {1365841,5},
1432       {2430194,2},
1433       {2147483647,4}
1434     }
1435     },
1436     { 64,8,{
1437       {0,3},
1438       {4,4},
1439       {17,5},
1440       {635,1},
1441       {202937,3},
1442       {308253,5},
1443       {1389874,4},
1444       {2147483647,2}
1445     }
1446     },
1447     { 128,8,{
1448       {0,3},
1449       {4,4},
1450       {16,5},
1451       {1238,1},
1452       {280097,3},
1453       {631434,5},
1454       {2605072,4},
1455       {2147483647,2}
1456     }
1457     },
1458     { 256,5,{
1459       {0,2},
1460       {4,4},
1461       {16,5},
1462       {2418,1},
1463       {2147483647,3}
1464     }
1465     },
1466     { 2147483647,6,{
1467       {0,2},
1468       {4,4},
1469       {16,5},
1470       {33182,1},
1471       {3763779,3},
1472       {2147483647,4}
1473     }
1474     }
1475   }
1476   }
1477 };
1478
1479 /* I_MPI_ADJUST_ALLGATHER
1480
1481 MPI_Allgather
1482
1483 1. Recursive doubling algorithm
1484 2. Bruck's algorithm
1485 3. Ring algorithm
1486 4. Topology aware Gatherv + Bcast algorithm
1487
1488 */
1489
1490 int (*intel_allgather_functions_table[])(void *sbuf, int scount,
1491                                               MPI_Datatype sdtype,
1492                                               void* rbuf, int rcount,
1493                                               MPI_Datatype rdtype,
1494                                               MPI_Comm  comm
1495                                                     ) ={
1496       Coll_allgather_rdb::allgather,
1497       Coll_allgather_bruck::allgather,
1498       Coll_allgather_ring::allgather,
1499       Coll_allgather_GB::allgather
1500 };
1501
1502 intel_tuning_table_element intel_allgather_table[] =
1503 {
1504   {1,{
1505     {4,11,
1506         {
1507         {1,4},
1508         {384,1},
1509         {1533,4},
1510         {3296,1},
1511         {10763,4},
1512         {31816,3},
1513         {193343,4},
1514         {405857,3},
1515         {597626,4},
1516         {1844323,3},
1517         {2147483647,4}
1518         }
1519     },
1520     {8,10,
1521         {
1522         {12,4},
1523         {46,1},
1524         {205,4},
1525         {3422,2},
1526         {4200,4},
1527         {8748,1},
1528         {24080,3},
1529         {33244,4},
1530         {371159,1},
1531         {2147483647,3}
1532         }
1533     },
1534     {16, 8,
1535         {
1536         {3,4},
1537         {53,1},
1538         {100,4},
1539         {170,1},
1540         {6077,4},
1541         {127644,1},
1542         {143741,4},
1543         {2147483647,3}
1544         }
1545     },
1546     {2147483647,10,
1547         {
1548         {184,1},
1549         {320,4},
1550         {759,1},
1551         {1219,4},
1552         {2633,1},
1553         {8259,4},
1554         {123678,1},
1555         {160801,4},
1556         {284341,1},
1557         {2147483647,4}
1558         }
1559     }
1560   }
1561   },
1562   {2,{
1563     { 8,6,{
1564       {490,1},
1565       {558,2},
1566       {2319,1},
1567       {46227,3},
1568       {2215101,1},
1569       {2147483647,3}
1570     }
1571     },
1572     { 16,4,{
1573       {1005,1},
1574       {1042,2},
1575       {2059,1},
1576       {2147483647,3}
1577     }
1578     },
1579     { 2147483647,2,{
1580       {2454,1},
1581       {2147483647,3}
1582    }
1583    }
1584   }
1585   },
1586   {4,{
1587     { 8,2,{
1588       {2861,1},
1589       {2147483647,3}
1590     }
1591     },
1592     { 2147483647,2,{
1593       {605,1},
1594       {2147483647,3}
1595     }
1596     }
1597   }
1598   },
1599   {8,{
1600     { 16,4,{
1601       {66,1},
1602       {213,4},
1603       {514,1},
1604       {2147483647,3}
1605     }
1606     },
1607     { 32,4,{
1608       {91,1},
1609       {213,4},
1610       {514,1},
1611       {2147483647,3}
1612     }
1613     },
1614     { 64,4,{
1615       {71,1},
1616       {213,4},
1617       {514,1},
1618       {2147483647,3}
1619     }
1620     },
1621     { 128,2,{
1622       {305,1},
1623       {2147483647,3}
1624     }
1625     },
1626     { 2147483647,2,{
1627       {213,1},
1628       {2147483647,3}
1629     }
1630     }
1631   }
1632   },
1633   {16,{
1634     { 8,2,{
1635       {2861,1},
1636       {2147483647,3}
1637     }
1638     },
1639     { 2147483647,2,{
1640       {605,1},
1641       {2147483647,3}
1642     }
1643     }
1644   }
1645   }
1646 };
1647
1648 /* I_MPI_ADJUST_ALLGATHERV
1649
1650 MPI_Allgatherv
1651
1652 1. Recursive doubling algorithm
1653 2. Bruck's algorithm
1654 3. Ring algorithm
1655 4. Topology aware Gatherv + Bcast algorithm
1656
1657 */
1658
1659 int (*intel_allgatherv_functions_table[])(void *sbuf, int scount,
1660                                                MPI_Datatype sdtype,
1661                                                void* rbuf, int *rcounts,
1662                                                int *rdispls,
1663                                                MPI_Datatype rdtype,
1664                                                MPI_Comm  comm
1665                                                     ) ={
1666       Coll_allgatherv_mpich_rdb::allgatherv,
1667       Coll_allgatherv_ompi_bruck::allgatherv,
1668       Coll_allgatherv_ring::allgatherv,
1669       Coll_allgatherv_GB::allgatherv
1670 };
1671
1672 intel_tuning_table_element intel_allgatherv_table[] =
1673 {
1674   {1,{
1675     {2,3,
1676         {
1677         {259668,3},
1678         {635750,4},
1679         {2147483647,3}
1680         }
1681     },
1682     {4,7,
1683         {
1684         {1,1},
1685         {5,4},
1686         {46,1},
1687         {2590,2},
1688         {1177259,3},
1689         {2767234,4},
1690         {2147483647,3}
1691         }
1692     },
1693     {8, 6,
1694         {
1695         {99,2},
1696         {143,1},
1697         {4646,2},
1698         {63522,3},
1699         {2187806,4},
1700         {2147483647,3}
1701         }
1702     },
1703     {2147483647,7,
1704         {
1705         {1,1},
1706         {5,4},
1707         {46,1},
1708         {2590,2},
1709         {1177259,3},
1710         {2767234,4},
1711         {2147483647,3}
1712         }
1713     }
1714   }
1715   },
1716   {2,{
1717     { 4,3,{
1718       {3147,1},
1719       {5622,2},
1720       {2147483647,3}
1721     }
1722     },
1723     { 8,3,{
1724       {975,1},
1725       {4158,2},
1726       {2147483647,3}
1727     }
1728     },
1729     { 16,2,{
1730       {2146,1},
1731       {2147483647,3}
1732     }
1733     },
1734     { 32,4,{
1735       {81,1},
1736       {414,2},
1737       {1190,1},
1738       {2147483647,3}
1739     }
1740     },
1741     { 2147483647,5,{
1742       {1,2},
1743       {3,1},
1744       {783,2},
1745       {1782,4},
1746       {2147483647,3}
1747     }
1748     }
1749   }
1750   },
1751   {4,{
1752     { 8,2,{
1753       {2554,1},
1754       {2147483647,3}
1755     }
1756     },
1757     { 16,4,{
1758       {272,1},
1759       {657,2},
1760       {2078,1},
1761       {2147483647,3}
1762     }
1763     },
1764     { 32,2,{
1765       {1081,1},
1766       {2147483647,3}
1767     }
1768     },
1769     { 64,2,{
1770       {547,1},
1771       {2147483647,3}
1772     }
1773     },
1774     { 2147483647,5,{
1775       {19,1},
1776       {239,2},
1777       {327,1},
1778       {821,4},
1779       {2147483647,3}
1780     }
1781     }
1782   }
1783   },
1784   {8,{
1785     { 16,3,{
1786       {55,1},
1787       {514,2},
1788       {2147483647,3}
1789     }
1790     },
1791     { 32,4,{
1792       {53,1},
1793       {167,4},
1794       {514,2},
1795       {2147483647,3}
1796     }
1797     },
1798     { 64,3,{
1799       {13,1},
1800       {319,4},
1801       {2147483647,3}
1802     }
1803     },
1804     { 128,7,{
1805       {2,1},
1806       {11,2},
1807       {48,1},
1808       {201,2},
1809       {304,1},
1810       {1048,4},
1811       {2147483647,3}
1812     }
1813     },
1814     { 2147483647,5,{
1815       {5,1},
1816       {115,4},
1817       {129,1},
1818       {451,4},
1819       {2147483647,3}
1820     }
1821     }
1822   }
1823   },
1824   {16,{
1825     { 8,2,{
1826       {2554,1},
1827       {2147483647,3}
1828     }
1829     },
1830     { 16,4,{
1831       {272,1},
1832       {657,2},
1833       {2078,1},
1834       {2147483647,3}
1835     }
1836     },
1837     { 32,2,{
1838       {1081,1},
1839       {2147483647,3}
1840     }
1841     },
1842     { 64,2,{
1843       {547,1},
1844       {2147483647,3}
1845     }
1846     },
1847     { 2147483647,5,{
1848       {19,1},
1849       {239,2},
1850       {327,1},
1851       {821,4},
1852       {2147483647,3}
1853     }
1854     }
1855   }
1856   }
1857 };
1858
1859
1860 /* I_MPI_ADJUST_GATHER
1861
1862 MPI_Gather
1863
1864 1. Binomial algorithm
1865 2. Topology aware binomial algorithm
1866 3. Shumilin's algorithm
1867
1868 */
1869
1870 int (*intel_gather_functions_table[])(void *sbuf, int scount,
1871                                            MPI_Datatype sdtype,
1872                                            void* rbuf, int rcount,
1873                                            MPI_Datatype rdtype,
1874                                            int root,
1875                                            MPI_Comm  comm
1876                                                     ) ={
1877       Coll_gather_ompi_binomial::gather,
1878       Coll_gather_ompi_binomial::gather,
1879       Coll_gather_mvapich2::gather
1880 };
1881
1882 intel_tuning_table_element intel_gather_table[] =
1883 {
1884   {1,{
1885     {8,3,
1886         {
1887         {17561,3},
1888         {44791,2},
1889         {2147483647,3}
1890         }
1891     },
1892     {16,7,
1893         {
1894         {16932,3},
1895         {84425,2},
1896         {158363,3},
1897         {702801,2},
1898         {1341444,3},
1899         {2413569,2},
1900         {2147483647,3}
1901         }
1902     },
1903     {2147483647,4,
1904         {
1905         {47187,3},
1906         {349696,2},
1907         {2147483647,3},
1908         {2147483647,1}
1909         }
1910     }
1911   }
1912   },
1913   {2,{
1914     {2147483647,1,{
1915       {2147483647,3}
1916     }
1917     }
1918   }
1919   },
1920   {4,{
1921     {2147483647,1,{
1922       {2147483647,3}
1923     }
1924     }
1925   }
1926   },
1927   {8,{
1928     { 16,1,{
1929       {2147483647,3}
1930     }
1931     },
1932     { 32,2,{
1933       {9,2},
1934       {2147483647,3}
1935     }
1936     },
1937     { 64,2,{
1938       {784,2},
1939       {2147483647,3}
1940     }
1941     },
1942     { 128,3,{
1943       {160,3},
1944       {655,2},
1945       {2147483647,3}
1946     }
1947     },
1948     { 2147483647,1,{
1949       {2147483647,3}
1950     }
1951     }
1952   }
1953   },
1954   {16,{
1955     {2147483647,1,{
1956       {2147483647,3}
1957     }
1958     }
1959   }
1960   }
1961 };
1962
1963
1964 /* I_MPI_ADJUST_SCATTER
1965
1966 MPI_Scatter
1967
1968 1. Binomial algorithm
1969 2. Topology aware binomial algorithm
1970 3. Shumilin's algorithm
1971
1972 */
1973
1974 int (*intel_scatter_functions_table[])(void *sbuf, int scount,
1975                                             MPI_Datatype sdtype,
1976                                             void* rbuf, int rcount,
1977                                             MPI_Datatype rdtype,
1978                                             int root, MPI_Comm  comm
1979                                                     ) ={
1980       Coll_scatter_ompi_binomial::scatter,
1981       Coll_scatter_ompi_binomial::scatter,
1982       Coll_scatter_mvapich2::scatter
1983 };
1984
1985 intel_tuning_table_element intel_scatter_table[] =
1986 {
1987   {1,{
1988     {2,2,
1989         {
1990         {16391,1},
1991         {2147483647,3}
1992         }
1993     },
1994     {4,6,
1995         {
1996         {16723,3},
1997         {153541,2},
1998         {425631,3},
1999         {794142,2},
2000         {1257027,3},
2001         {2147483647,2}
2002         }
2003     },
2004     {8,7,
2005         {
2006         {2633,3},
2007         {6144,2},
2008         {14043,3},
2009         {24576,2},
2010         {107995,3},
2011         {1752729,2},
2012         {2147483647,3}
2013         }
2014     },
2015     {16,7,
2016         {
2017         {2043,3},
2018         {2252,2},
2019         {17749,3},
2020         {106020,2},
2021         {628654,3},
2022         {3751354,2},
2023         {2147483647,3}
2024         }
2025     },
2026     {2147483647,4,
2027         {
2028         {65907,3},
2029         {245132,2},
2030         {1042439,3},
2031         {2147483647,2},
2032         {2147483647,1}
2033         }
2034     }
2035   }
2036   },
2037   {2,{
2038      {2147483647,1,{
2039        {2147483647,3}
2040      }
2041      }
2042   }
2043   },
2044   {4,{
2045     { 8,1,{
2046       {2147483647,3}
2047     }
2048     },
2049     { 16,2,{
2050       {140,3},
2051       {1302,1},
2052       {2147483647,3}
2053     }
2054     },
2055     { 32,2,{
2056       {159,3},
2057       {486,1},
2058       {2147483647,3}
2059     }
2060     },
2061     { 64,2,{
2062       {149,1},
2063       {2147483647,3}
2064     }
2065     },
2066     { 2147483647,2,{
2067       {139,1},
2068       {2147483647,3}
2069     }
2070     }
2071   }
2072   },
2073   {8,{
2074     { 16,4,{
2075       {587,1},
2076       {1370,2},
2077       {2102,1},
2078       {2147483647,3}
2079     }
2080     },
2081     { 32,3,{
2082       {1038,1},
2083       {2065,2},
2084       {2147483647,3}
2085     }
2086     },
2087     { 64,3,{
2088       {515,1},
2089       {2069,2},
2090       {2147483647,3}
2091     }
2092     },
2093     { 128,3,{
2094       {284,1},
2095       {796,2},
2096       {2147483647,3}
2097     }
2098     },
2099     { 2147483647,2,{
2100       {139,1},
2101       {2147483647,3}
2102     }
2103     }
2104   }
2105   },
2106   {16,{
2107     { 8,1,{
2108       {2147483647,3}
2109     }
2110     },
2111     { 16,3,{
2112       {140,3},
2113       {1302,1},
2114       {2147483647,3}
2115     }
2116     },
2117     { 32,3,{
2118       {159,3},
2119       {486,1},
2120       {2147483647,3}
2121     }
2122     },
2123     { 64,2,{
2124       {149,1},
2125       {2147483647,3}
2126     }
2127     },
2128     { 2147483647,2,{
2129       {139,1},
2130       {2147483647,3}
2131     }
2132     }
2133   }
2134   }
2135 };
2136
2137
2138
2139 /* I_MPI_ADJUST_ALLTOALLV
2140
2141 MPI_Alltoallv
2142
2143 1. Isend/Irecv + waitall algorithm
2144 2. Plum's algorithm
2145
2146 */
2147
2148 int (*intel_alltoallv_functions_table[])(void *sbuf, int *scounts, int *sdisps,
2149                                               MPI_Datatype sdtype,
2150                                               void *rbuf, int *rcounts, int *rdisps,
2151                                               MPI_Datatype rdtype,
2152                                               MPI_Comm  comm
2153                                                     ) ={
2154       Coll_alltoallv_ompi_basic_linear::alltoallv,
2155       Coll_alltoallv_bruck::alltoallv
2156 };
2157
2158 intel_tuning_table_element intel_alltoallv_table[] =
2159 {
2160   {1,{
2161     {2147483647,1,
2162         {
2163         {2147483647,1}
2164         }
2165     }
2166   }
2167   },
2168   {2,{
2169     {2147483647,1,
2170         {
2171         {2147483647,1}
2172         }
2173     }
2174   }
2175   },
2176   {4,{
2177     { 8,1,{
2178       {2147483647,1}//weirdly, intel reports the use of algo 0 here
2179     }
2180     },
2181     { 2147483647,2,{
2182       {4,1},//0 again
2183       {2147483647,2}
2184     }
2185     }
2186   }
2187   },
2188   {8,{
2189     { 16,1,{
2190       {2147483647,1}
2191     }
2192     },
2193     { 2147483647,2,{
2194       {0,1},//weird again, only for 0-sized messages
2195       {2147483647,2}
2196     }
2197     }
2198   }
2199   },
2200   {16,{
2201     { 8,1,{
2202       {2147483647,1}//0
2203     }
2204     },
2205     { 2147483647,2,{
2206       {4,1},//0
2207       {2147483647,2}
2208     }
2209     }
2210   }
2211   }
2212 };
2213
2214
2215 //These are collected from table 3.5-2 of the Intel MPI Reference Manual
2216
2217
2218 #define SIZECOMP_reduce_scatter\
2219     int total_message_size = 0;\
2220     for (i = 0; i < comm_size; i++) { \
2221         total_message_size += rcounts[i];\
2222     }\
2223     size_t block_dsize = total_message_size*dtype->size();\
2224
2225 #define SIZECOMP_allreduce\
2226   size_t block_dsize =rcount * dtype->size();
2227
2228 #define SIZECOMP_alltoall\
2229   size_t block_dsize =send_count * send_type->size();
2230
2231 #define SIZECOMP_bcast\
2232   size_t block_dsize =count * datatype->size();
2233
2234 #define SIZECOMP_reduce\
2235   size_t block_dsize =count * datatype->size();
2236
2237 #define SIZECOMP_barrier\
2238   size_t block_dsize = 1;
2239
2240 #define SIZECOMP_allgather\
2241   size_t block_dsize =recv_count * recv_type->size();
2242
2243 #define SIZECOMP_allgatherv\
2244     int total_message_size = 0;\
2245     for (i = 0; i < comm_size; i++) { \
2246         total_message_size += recv_count[i];\
2247     }\
2248     size_t block_dsize = total_message_size*recv_type->size();
2249
2250 #define SIZECOMP_gather\
2251   int rank = comm->rank();\
2252   size_t block_dsize = (send_buff == MPI_IN_PLACE || rank ==root) ?\
2253                 recv_count * recv_type->size() :\
2254                 send_count * send_type->size();
2255
2256 #define SIZECOMP_scatter\
2257   int rank = comm->rank();\
2258   size_t block_dsize = (sendbuf == MPI_IN_PLACE || rank !=root ) ?\
2259                 recvcount * recvtype->size() :\
2260                 sendcount * sendtype->size();
2261
2262 #define SIZECOMP_alltoallv\
2263   size_t block_dsize = 1;
2264
2265 #define IMPI_COLL_SELECT(cat, ret, args, args2)\
2266 ret Coll_ ## cat ## _impi:: cat (COLL_UNPAREN args)\
2267 {\
2268     int comm_size = comm->size();\
2269     int i =0;\
2270     SIZECOMP_ ## cat\
2271     i=0;\
2272     int j =0, k=0;\
2273     if(comm->get_leaders_comm()==MPI_COMM_NULL){\
2274       comm->init_smp();\
2275     }\
2276     int local_size=1;\
2277     if (comm->is_uniform()) {\
2278         local_size = comm->get_intra_comm()->size();\
2279     }\
2280     while(i < INTEL_MAX_NB_PPN &&\
2281     local_size!=intel_ ## cat ## _table[i].ppn)\
2282       i++;\
2283     if(i==INTEL_MAX_NB_PPN) i=0;\
2284     while(comm_size>intel_ ## cat ## _table[i].elems[j].max_num_proc\
2285         && j < INTEL_MAX_NB_THRESHOLDS)\
2286       j++;\
2287     while(block_dsize >=intel_ ## cat ## _table[i].elems[j].elems[k].max_size\
2288          && k< intel_ ## cat ## _table[i].elems[j].num_elems)\
2289       k++;\
2290     return (intel_ ## cat ## _functions_table[intel_ ## cat ## _table[i].elems[j].elems[k].algo-1]\
2291     args2);\
2292 }
2293
2294
2295 COLL_APPLY(IMPI_COLL_SELECT, COLL_ALLGATHERV_SIG, (send_buff, send_count, send_type, recv_buff, recv_count, recv_disps, recv_type, comm));
2296 COLL_APPLY(IMPI_COLL_SELECT, COLL_ALLREDUCE_SIG, (sbuf, rbuf, rcount, dtype, op, comm));
2297 COLL_APPLY(IMPI_COLL_SELECT, COLL_GATHER_SIG, (send_buff, send_count, send_type, recv_buff, recv_count, recv_type, root, comm));
2298 COLL_APPLY(IMPI_COLL_SELECT, COLL_ALLGATHER_SIG, (send_buff,send_count,send_type,recv_buff,recv_count,recv_type,comm));
2299 COLL_APPLY(IMPI_COLL_SELECT, COLL_ALLTOALL_SIG,(send_buff, send_count, send_type, recv_buff, recv_count, recv_type,comm));
2300 COLL_APPLY(IMPI_COLL_SELECT, COLL_ALLTOALLV_SIG, (send_buff, send_counts, send_disps, send_type, recv_buff, recv_counts, recv_disps, recv_type, comm));
2301 COLL_APPLY(IMPI_COLL_SELECT, COLL_BCAST_SIG , (buf, count, datatype, root, comm));
2302 COLL_APPLY(IMPI_COLL_SELECT, COLL_REDUCE_SIG,(buf,rbuf, count, datatype, op, root, comm));
2303 COLL_APPLY(IMPI_COLL_SELECT, COLL_REDUCE_SCATTER_SIG ,(sbuf,rbuf, rcounts,dtype,op,comm));
2304 COLL_APPLY(IMPI_COLL_SELECT, COLL_SCATTER_SIG ,(sendbuf, sendcount, sendtype,recvbuf, recvcount, recvtype,root, comm));
2305 COLL_APPLY(IMPI_COLL_SELECT, COLL_BARRIER_SIG,(comm));
2306
2307 }
2308 }