Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
avoid word being recognized as special by doxygen
[simgrid.git] / doc / doxygen / module-smpi.doc
1 /** @addtogroup SMPI_API
2
3 This programming environment enables the study of MPI application by
4 emulating them on top of the SimGrid simulator. This is particularly
5 interesting to study existing MPI applications within the comfort of
6 the simulator. The motivation for this work is detailed in the
7 reference article (available at http://hal.inria.fr/inria-00527150).
8
9 Our goal is to enable the study of *unmodified* MPI applications,
10 although we are not quite there yet (see @ref SMPI_what). In
11 addition, you can modify your code to speed up your studies or
12 otherwise increase their scalability (see @ref SMPI_adapting).
13
14 \section SMPI_who Who should use SMPI (and who shouldn't)
15
16 SMPI is now considered as stable and you can use it in production. You
17 may probably want to read the scientific publications that detail the
18 models used and their limits, but this should not be absolutely
19 necessary. If you already fluently write and use MPI applications,
20 SMPI should sound very familiar to you. Use smpicc instead of mpicc,
21 and smpirun instead of mpirun (see below for more details).
22
23 Of course, if you don't know what MPI is, the documentation of SMPI
24 will seem a bit terse to you. You should pick up a good MPI tutorial
25 on the Internet (or a course in your favorite university) and come
26 back to SMPI once you know a bit more about MPI. Alternatively, you
27 may want to turn to the other SimGrid interfaces such as the 
28 \ref MSG_API environment, or the \ref SD_API one.
29
30 \section SMPI_what What can run within SMPI?
31
32 You can run unmodified MPI applications (both C and Fortran) within
33 SMPI, provided that (1) you only use MPI calls that we implemented in
34 MPI and (2) you don't use any globals in your application.
35
36 \subsection SMPI_what_coverage MPI coverage of SMPI
37
38 Our coverage of the interface is very decent, but still incomplete;
39 Given the size of the MPI standard, it may well be that we never
40 implement absolutely all existing primitives. One sided communications
41 and I/O primitives are not targeted for now. Our current state is
42 still very decent: we pass most of the MPICH coverage tests.
43
44 The full list of not yet implemented functions is documented in the
45 file <tt>include/smpi/smpi.h</tt> of the archive, between two lines
46 containing the <tt>FIXME</tt> marker. If you really need a missing
47 feature, please get in touch with us: we can guide you though the
48 SimGrid code to help you implementing it, and we'd glad to integrate
49 it in the main project afterward if you contribute them back.
50
51 \subsection SMPI_what_globals Global variables
52
53 Concerning the globals, the problem comes from the fact that usually,
54 MPI processes run as real UNIX processes while they are all folded
55 into threads of a unique system process in SMPI. Global variables are
56 usually private to each MPI process while they become shared between
57 the processes in SMPI. This point is rather problematic, and currently
58 forces to modify your application to privatize the global variables.
59
60 We tried several techniques to work this around. We used to have a
61 script that privatized automatically the globals through static
62 analysis of the source code, but it was not robust enough to be used
63 in production. This issue, as well as several potential solutions, is
64 discussed in this article: "Automatic Handling of Global Variables for
65 Multi-threaded MPI Programs",
66 available at http://charm.cs.illinois.edu/newPapers/11-23/paper.pdf
67 (note that this article does not deal with SMPI but with a concurrent
68 solution called AMPI that suffers of the same issue). 
69
70 A method using dynamic switching of the .data and .bss segments of an
71 ELF executable has been introduced in SimGrid 3.11. By using the <tt>smpi/
72 privatize_global_variableles</tt> option to yes, SMPI will duplicate
73 the segments containing the global variables and when needed, will map 
74 the right one in memory. This needs ELF executables and mmap on the system
75 (Linux and recent BSDs should be compatible). %As no copy is involved, 
76 performance should not be altered (but memory occupation will be higher).
77
78 This solution actually works really good for a good number of MPI 
79 applications. Its main limitation is that if the application loads dynamic 
80 libraries, their global variables won't be privatized. This can be avoided 
81 by linking statically with these libraries (but NOT with libsimgrid, as we 
82 need SimGrid's own global varibles).
83
84
85 \section SMPI_compiling Compiling your code
86
87 This is very simply done with the <tt>smpicc</tt> script. If you
88 already compiled any MPI code before, you already know how to use it.
89 If not, you should try to get your MPI code running on top of MPI
90 before giving SMPI a spin. Actually, that's very simple even if it's
91 the first time you use MPI code: just use smpicc as a compiler (in
92 replacement of gcc or your usual compiler), and you're set.
93
94 \section SMPI_executing Executing your code on top of the simulator
95
96 This is done though the <tt>smpirun</tt> script as follows.
97 <tt>my_hostfile.txt</tt> is a classical MPI hostfile (that is, this
98 file lists the machines on which the processes must be dispatched, one
99 per line)  <tt>my_platform.xml</tt> is a classical SimGrid platform
100 file. Of course, the hosts of the hostfile must exist in the provided
101 platform. <tt>./program</tt> is the MPI program that you want to
102 simulate (must be compiled by <tt>smpicc</tt>) while <tt>-arg</tt> is
103 a command-line parameter passed to this program.
104
105 \verbatim
106 smpirun -hostfile my_hostfile.txt -platform my_platform.xml ./program -arg
107 \endverbatim
108
109 smpirun accepts other parameters, such as <tt>-np</tt> if you don't
110 want to use all the hosts defined in the hostfile, <tt>-map</tt> to
111 display on which host each rank gets mapped of <tt>-trace</tt> to
112 activate the tracing during the simulation. You can get the full list
113 by running
114 \verbatim
115 smpirun -help
116 \endverbatim
117
118 \section SMPI_adapting Adapting your MPI code to the use of SMPI
119
120 As detailed in the reference article (available at
121 http://hal.inria.fr/inria-00527150), you may want to adapt your code
122 to improve the simulation performance. But these tricks may seriously
123 hinder the result quality (or even prevent the app to run) if used
124 wrongly. We assume that if you want to simulate an HPC application,
125 you know what you are doing. Don't prove us wrong!
126
127 \section SMPI_adapting_size Reducing your memory footprint
128
129 If you get short on memory (the whole app is executed on a single node when
130 simulated), you should have a look at the SMPI_SHARED_MALLOC and
131 SMPI_SHARED_FREE macros. It allows to share memory areas between processes: The
132 purpose of these macro is that the same line malloc on each process will point
133 to the exact same memory area. So if you have a malloc of 2M and you have 16
134 processes, this macro will change your memory consumption from 2M*16 to 2M
135 only. Only one block for all processes.
136
137 If your program is ok with a block containing garbage value because all
138 processes write and read to the same place without any kind of coordination,
139 then this macro can dramatically shrink your memory consumption. For example,
140 that will be very beneficial to a matrix multiplication code, as all blocks will
141 be stored on the same area. Of course, the resulting computations will useless,
142 but you can still study the application behavior this way. 
143
144 Naturally, this won't work if your code is data-dependent. For example, a Jacobi
145 iterative computation depends on the result computed by the code to detect
146 convergence conditions, so turning them into garbage by sharing the same memory
147 area between processes does not seem very wise. You cannot use the
148 SMPI_SHARED_MALLOC macro in this case, sorry.
149
150 This feature is demoed by the example file
151 <tt>examples/smpi/NAS/DT-folding/dt.c</tt>
152
153 \section SMPI_adapting_speed Toward faster simulations
154
155 If your application is too slow, try using SMPI_SAMPLE_LOCAL,
156 SMPI_SAMPLE_GLOBAL and friends to indicate which computation loops can
157 be sampled. Some of the loop iterations will be executed to measure
158 their duration, and this duration will be used for the subsequent
159 iterations. These samples are done per processor with
160 SMPI_SAMPLE_LOCAL, and shared between all processors with
161 SMPI_SAMPLE_GLOBAL. Of course, none of this will work if the execution
162 time of your loop iteration are not stable.
163
164 This feature is demoed by the example file 
165 <tt>examples/smpi/NAS/EP-sampling/ep.c</tt>
166
167
168 \section SMPI_collective_algorithms Simulating collective operations
169
170 MPI collective operations can be implemented very differently from one library 
171 to another. Actually, all existing libraries implement several algorithms 
172 for each collective operation, and by default select at runtime which one 
173 should be used for the current operation, depending on the sizes sent, the number
174  of nodes, the communicator, or the communication library being used. These 
175 decisions are based on empirical results and theoretical complexity estimation, 
176 but they can sometimes be suboptimal. Manual selection is possible in these cases, 
177 to allow the user to tune the library and use the better collective if the 
178 default one is not good enough.
179
180 SMPI tries to apply the same logic, regrouping algorithms from OpenMPI, MPICH 
181 libraries, StarMPI (<a href="http://star-mpi.sourceforge.net/">STAR-MPI</a>), and MVAPICH2 libraries.
182 This collection of more than 115 algorithms allows a simple and effective
183  comparison of their behavior and performance, making SMPI a tool of choice for the
184 development of such algorithms.
185
186 \subsection Tracing_internals Tracing of internal communications
187
188 For each collective, default tracing only outputs global data. 
189 Internal communication operations are not traced to avoid outputting too much data
190 to the trace. To debug and compare algorithm, this can be changed with the item 
191 \b tracing/smpi/internals , which has 0 for default value.
192 Here are examples of two alltoall collective algorithms runs on 16 nodes, 
193 the first one with a ring algorithm, the second with a pairwise one :
194
195 \htmlonly
196 <a href="smpi_simgrid_alltoall_ring_16.png" border=0><img src="smpi_simgrid_alltoall_ring_16.png" width="30%" border=0 align="center"></a>
197 <a href="smpi_simgrid_alltoall_pair_16.png" border=0><img src="smpi_simgrid_alltoall_pair_16.png" width="30%" border=0 align="center"></a>
198 <br/>
199 \endhtmlonly
200
201 \subsection Selectors
202
203 The default selection logic implemented by default in OpenMPI (version 1.7) 
204 and MPICH (version 3.0.4) has been replicated and can be used by setting the
205 \b smpi/coll_selector item to either ompi or mpich. A selector based on the selection logic of MVAPICH2 (version 1.9) tuned on the Stampede cluster as also been implemented, as well as a preliminary version of an Intel MPI selector (version 4.1.3, also tuned for the Stampede cluster). Due the closed source nature of Intel MPI, some of the algorithms described in the documentation are not available, and are replaced by mvapich ones.
206
207 Values for option \b smpi/coll_selector are :
208  - ompi
209  - mpich
210  - mvapich2
211  - impi
212  - default
213
214 The code and details for each 
215 selector can be found in the <tt>src/smpi/colls/smpi_(openmpi/mpich/mvapich2/impi)_selector.c</tt> file.
216 As this is still in development, we do not insure that all algorithms are correctly
217  replicated and that they will behave exactly as the real ones. If you notice a difference,
218 please contact <a href="http://lists.gforge.inria.fr/mailman/listinfo/simgrid-devel">SimGrid developers mailing list</a>
219
220 The default selector uses the legacy algorithms used in versions of SimGrid
221  previous to the 3.10. they should not be used to perform performance study and 
222 may be removed in the future, a different selector being used by default.
223
224 \subsection algos Available algorithms
225
226 For each one of the listed algorithms, several versions are available,
227  either coming from STAR-MPI, MPICH or OpenMPI implementations. Details can be
228  found in the code or in <a href="http://www.cs.arizona.edu/~dkl/research/papers/ics06.pdf">STAR-MPI</a> for STAR-MPI algorithms.
229
230 Each collective can be selected using the corresponding configuration item. For example, to use the pairwise alltoall algorithm, one should add \b --cfg=smpi/alltoall:pair to the line. This will override the selector (for this algorithm only) if provided, allowing better flexibility.
231
232 Warning: Some collective may require specific conditions to be executed correctly (for instance having a communicator with a power of two number of nodes only), which are currently not enforced by Simgrid. Some crashes can be expected while trying these algorithms with unusual sizes/parameters
233
234 \subsubsection MPI_Alltoall
235
236 Most of these are best described in <a href="http://www.cs.arizona.edu/~dkl/research/papers/ics06.pdf">STAR-MPI</a>
237
238  - default : naive one, by default
239  - ompi : use openmpi selector for the alltoall operations
240  - mpich : use mpich selector for the alltoall operations
241  - mvapich2 : use mvapich2 selector for the alltoall operations
242  - impi : use intel mpi selector for the alltoall operations
243  - automatic (experimental) : use an automatic self-benchmarking algorithm 
244  - 2dmesh : organizes the nodes as a two dimensional mesh, and perform allgather 
245 along the dimensions
246  - 3dmesh : adds a third dimension to the previous algorithm
247  - rdb : recursive doubling : extends the mesh to a nth dimension, each one 
248 containing two nodes
249  - pair : pairwise exchange, only works for power of 2 procs, size-1 steps,
250 each process sends and receives from the same process at each step
251  - pair_light_barrier : same, with small barriers between steps to avoid contention
252  - pair_mpi_barrier : same, with MPI_Barrier used
253  - pair_one_barrier : only one barrier at the beginning
254  - ring : size-1 steps, at each step a process send to process (n+i)%size, and receives from (n-i)%size
255  - ring_light_barrier : same, with small barriers between some phases to avoid contention
256  - ring_mpi_barrier : same, with MPI_Barrier used
257  - ring_one_barrier : only one barrier at the beginning
258  - basic_linear : posts all receives and all sends,
259 starts the communications, and waits for all communication to finish
260  - mvapich2_scatter_dest : isend/irecv with scattered destinations, posting only a few messages at the same time
261
262 \subsubsection MPI_Alltoallv
263
264  - default : naive one, by default
265  - ompi : use openmpi selector for the alltoallv operations
266  - mpich : use mpich selector for the alltoallv operations
267  - mvapich2 : use mvapich2 selector for the alltoallv operations
268  - impi : use intel mpi selector for the alltoallv operations
269  - automatic (experimental) : use an automatic self-benchmarking algorithm 
270  - bruck : same as alltoall
271  - pair : same as alltoall
272  - pair_light_barrier : same as alltoall
273  - pair_mpi_barrier : same as alltoall
274  - pair_one_barrier : same as alltoall
275  - ring : same as alltoall
276  - ring_light_barrier : same as alltoall
277  - ring_mpi_barrier : same as alltoall
278  - ring_one_barrier : same as alltoall
279  - ompi_basic_linear : same as alltoall
280
281
282 \subsubsection MPI_Gather
283
284  - default : naive one, by default
285  - ompi : use openmpi selector for the gather operations
286  - mpich : use mpich selector for the gather operations
287  - mvapich2 : use mvapich2 selector for the gather operations
288  - impi : use intel mpi selector for the gather operations
289  - automatic (experimental) : use an automatic self-benchmarking algorithm 
290 which will iterate over all implemented versions and output the best
291  - ompi_basic_linear : basic linear algorithm from openmpi, each process sends to the root
292  - ompi_binomial : binomial tree algorithm
293  - ompi_linear_sync : same as basic linear, but with a synchronization at the
294  beginning and message cut into two segments.
295  - mvapich2_two_level : SMP-aware version from MVAPICH. Gather first intra-node (defaults to mpich's gather), and then exchange with only one process/node. Use mvapich2 selector to change these to tuned algorithms for Stampede cluster.
296
297 \subsubsection MPI_Barrier
298  - default : naive one, by default
299  - ompi : use openmpi selector for the barrier operations
300  - mpich : use mpich selector for the barrier operations
301  - mvapich2 : use mvapich2 selector for the barrier operations
302  - impi : use intel mpi selector for the barrier operations
303  - automatic (experimental) : use an automatic self-benchmarking algorithm 
304  - ompi_basic_linear : all processes send to root
305  - ompi_two_procs : special case for two processes
306  - ompi_bruck : nsteps = sqrt(size), at each step, exchange data with rank-2^k and rank+2^k
307  - ompi_recursivedoubling : recursive doubling algorithm
308  - ompi_tree : recursive doubling type algorithm, with tree structure
309  - ompi_doublering : double ring algorithm
310  - mvapich2_pair : pairwise algorithm
311
312
313 \subsubsection MPI_Scatter
314  - default : naive one, by default
315  - ompi : use openmpi selector for the scatter operations
316  - mpich : use mpich selector for the scatter operations
317  - mvapich2 : use mvapich2 selector for the scatter operations
318  - impi : use intel mpi selector for the scatter operations
319  - automatic (experimental) : use an automatic self-benchmarking algorithm 
320  - ompi_basic_linear : basic linear scatter 
321  - ompi_binomial : binomial tree scatter
322  - mvapich2_two_level_direct : SMP aware algorithm, with an intra-node stage (default set to mpich selector), and then a basic linear inter node stage. Use mvapich2 selector to change these to tuned algorithms for Stampede cluster. 
323  - mvapich2_two_level_binomial : SMP aware algorithm, with an intra-node stage (default set to mpich selector), and then a binomial phase. Use mvapich2 selector to change these to tuned algorithms for Stampede cluster.
324
325
326
327 \subsubsection MPI_Reduce
328  - default : naive one, by default
329  - ompi : use openmpi selector for the reduce operations
330  - mpich : use mpich selector for the reduce operations
331  - mvapich2 : use mvapich2 selector for the reduce operations
332  - impi : use intel mpi selector for the reduce operations
333  - automatic (experimental) : use an automatic self-benchmarking algorithm 
334  - arrival_pattern_aware : root exchanges with the first process to arrive
335  - binomial : uses a binomial tree
336  - flat_tree : uses a flat tree
337  - NTSL : Non-topology-specific pipelined linear-bcast function 
338    0->1, 1->2 ,2->3, ....., ->last node : in a pipeline fashion, with segments
339  of 8192 bytes
340  - scatter_gather : scatter then gather
341  - ompi_chain : openmpi reduce algorithms are built on the same basis, but the
342  topology is generated differently for each flavor
343 chain = chain with spacing of size/2, and segment size of 64KB 
344  - ompi_pipeline : same with pipeline (chain with spacing of 1), segment size 
345 depends on the communicator size and the message size
346  - ompi_binary : same with binary tree, segment size of 32KB
347  - ompi_in_order_binary : same with binary tree, enforcing order on the 
348 operations
349  - ompi_binomial : same with binomial algo (redundant with default binomial 
350 one in most cases)
351  - ompi_basic_linear : basic algorithm, each process sends to root
352  - mvapich2_knomial : k-nomial algorithm. Default factor is 4 (mvapich2 selector adapts it through tuning)
353  - mvapich2_two_level : SMP-aware reduce, with default set to mpich both for intra and inter communicators. Use mvapich2 selector to change these to tuned algorithms for Stampede cluster.
354  - rab : <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a>'s reduce algorithm 
355
356 \subsubsection MPI_Allreduce
357  - default : naive one, by default
358  - ompi : use openmpi selector for the allreduce operations
359  - mpich : use mpich selector for the allreduce operations
360  - mvapich2 : use mvapich2 selector for the allreduce operations
361  - impi : use intel mpi selector for the allreduce operations
362  - automatic (experimental) : use an automatic self-benchmarking algorithm 
363  - lr : logical ring reduce-scatter then logical ring allgather
364  - rab1 : variations of the  <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> algorithm : reduce_scatter then allgather
365  - rab2 : variations of the  <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> algorithm : alltoall then allgather
366  - rab_rsag : variation of the  <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> algorithm : recursive doubling 
367 reduce_scatter then recursive doubling allgather 
368  - rdb : recursive doubling
369  - smp_binomial : binomial tree with smp : binomial intra 
370 SMP reduce, inter reduce, inter broadcast then intra broadcast
371  - smp_binomial_pipeline : same with segment size = 4096 bytes
372  - smp_rdb : intra : binomial allreduce, inter : Recursive 
373 doubling allreduce, intra : binomial broadcast
374  - smp_rsag : intra : binomial allreduce, inter : reduce-scatter, 
375 inter:allgather, intra : binomial broadcast
376  - smp_rsag_lr : intra : binomial allreduce, inter : logical ring 
377 reduce-scatter, logical ring inter:allgather, intra : binomial broadcast
378  - smp_rsag_rab : intra : binomial allreduce, inter : rab
379 reduce-scatter, rab inter:allgather, intra : binomial broadcast
380  - redbcast : reduce then broadcast, using default or tuned algorithms if specified
381  - ompi_ring_segmented : ring algorithm used by OpenMPI
382  - mvapich2_rs : rdb for small messages, reduce-scatter then allgather else
383  - mvapich2_two_level : SMP-aware algorithm, with mpich as intra algoritm, and rdb as inter (Change this behavior by using mvapich2 selector to use tuned values)
384  - rab : default <a href="https://fs.hlrs.de/projects/par/mpi//myreduce.html">Rabenseifner</a> implementation
385
386 \subsubsection MPI_Reduce_scatter
387  - default : naive one, by default
388  - ompi : use openmpi selector for the reduce_scatter operations
389  - mpich : use mpich selector for the reduce_scatter operations
390  - mvapich2 : use mvapich2 selector for the reduce_scatter operations
391  - impi : use intel mpi selector for the reduce_scatter operations
392  - automatic (experimental) : use an automatic self-benchmarking algorithm 
393  - ompi_basic_recursivehalving : recursive halving version from OpenMPI
394  - ompi_ring : ring version from OpenMPI
395  - mpich_pair : pairwise exchange version from MPICH
396  - mpich_rdb : recursive doubling version from MPICH
397  - mpich_noncomm : only works for power of 2 procs, recursive doubling for noncommutative ops
398
399
400 \subsubsection MPI_Allgather
401
402  - default : naive one, by default
403  - ompi : use openmpi selector for the allgather operations
404  - mpich : use mpich selector for the allgather operations
405  - mvapich2 : use mvapich2 selector for the allgather operations
406  - impi : use intel mpi selector for the allgather operations
407  - automatic (experimental) : use an automatic self-benchmarking algorithm 
408  - 2dmesh : see alltoall
409  - 3dmesh : see alltoall
410  - bruck : Described by Bruck et.al. in <a href="http://ieeexplore.ieee.org/xpl/articleDetails.jsp?arnumber=642949">
411 Efficient algorithms for all-to-all communications in multiport message-passing systems</a> 
412  - GB : Gather - Broadcast (uses tuned version if specified)
413  - loosely_lr : Logical Ring with grouping by core (hardcoded, default 
414 processes/node: 4)
415  - NTSLR : Non Topology Specific Logical Ring
416  - NTSLR_NB : Non Topology Specific Logical Ring, Non Blocking operations
417  - pair : see alltoall
418  - rdb : see alltoall
419  - rhv : only power of 2 number of processes
420  - ring : see alltoall
421  - SMP_NTS : gather to root of each SMP, then every root of each SMP node 
422 post INTER-SMP Sendrecv, then do INTRA-SMP Bcast for each receiving message, 
423 using logical ring algorithm (hardcoded, default processes/SMP: 8)
424  - smp_simple : gather to root of each SMP, then every root of each SMP node 
425 post INTER-SMP Sendrecv, then do INTRA-SMP Bcast for each receiving message, 
426 using simple algorithm (hardcoded, default processes/SMP: 8)
427  - spreading_simple : from node i, order of communications is i -> i + 1, i ->
428  i + 2, ..., i -> (i + p -1) % P
429  - ompi_neighborexchange : Neighbor Exchange algorithm for allgather. 
430 Described by Chen et.al. in  <a href="http://ieeexplore.ieee.org/xpl/articleDetails.jsp?tp=&arnumber=1592302">Performance Evaluation of Allgather Algorithms on Terascale Linux Cluster with Fast Ethernet</a>
431  - mvapich2_smp : SMP aware algorithm, performing intra-node gather, inter-node allgather with one process/node, and bcast intra-node
432
433
434 \subsubsection MPI_Allgatherv
435  - default : naive one, by default
436  - ompi : use openmpi selector for the allgatherv operations
437  - mpich : use mpich selector for the allgatherv operations
438  - mvapich2 : use mvapich2 selector for the allgatherv operations
439  - impi : use intel mpi selector for the allgatherv operations
440  - automatic (experimental) : use an automatic self-benchmarking algorithm 
441  - GB : Gatherv - Broadcast (uses tuned version if specified, but only for 
442 Bcast, gatherv is not tuned)
443  - pair : see alltoall
444  - ring : see alltoall
445  - ompi_neighborexchange : see allgather
446  - ompi_bruck : see allgather
447  - mpich_rdb : recursive doubling algorithm from MPICH
448  - mpich_ring : ring algorithm from MPICh - performs differently from the 
449 one from STAR-MPI
450
451 \subsubsection MPI_Bcast
452  - default : naive one, by default
453  - ompi : use openmpi selector for the bcast operations
454  - mpich : use mpich selector for the bcast operations
455  - mvapich2 : use mvapich2 selector for the bcast operations
456  - impi : use intel mpi selector for the bcast operations
457  - automatic (experimental) : use an automatic self-benchmarking algorithm 
458  - arrival_pattern_aware : root exchanges with the first process to arrive
459  - arrival_pattern_aware_wait : same with slight variation
460  - binomial_tree : binomial tree exchange
461  - flattree : flat tree exchange
462  - flattree_pipeline : flat tree exchange, message split into 8192 bytes pieces
463  - NTSB : Non-topology-specific pipelined binary tree with 8192 bytes pieces
464  - NTSL : Non-topology-specific pipelined linear with 8192 bytes pieces
465  - NTSL_Isend : Non-topology-specific pipelined linear with 8192 bytes pieces, asynchronous communications
466  - scatter_LR_allgather : scatter followed by logical ring allgather
467  - scatter_rdb_allgather : scatter followed by recursive doubling allgather
468  - arrival_scatter : arrival pattern aware scatter-allgather
469  - SMP_binary : binary tree algorithm with 8 cores/SMP
470  - SMP_binomial : binomial tree algorithm with 8 cores/SMP
471  - SMP_linear : linear algorithm with 8 cores/SMP
472  - ompi_split_bintree : binary tree algorithm from OpenMPI, with message split in 8192 bytes pieces
473  - ompi_pipeline : pipeline algorithm from OpenMPI, with message split in 128KB pieces
474  - mvapich2_inter_node : Inter node default mvapich worker 
475  - mvapich2_intra_node : Intra node default mvapich worker
476  - mvapich2_knomial_intra_node :  k-nomial intra node default mvapich worker. default factor is 4.
477
478 \subsection auto Automatic evaluation 
479
480 (Warning : This is experimental and may be removed or crash easily)
481
482 An automatic version is available for each collective (or even as a selector). This specific 
483 version will loop over all other implemented algorithm for this particular collective, and apply 
484 them while benchmarking the time taken for each process. It will then output the quickest for 
485 each process, and the global quickest. This is still unstable, and a few algorithms which need 
486 specific number of nodes may crash.
487
488
489 \subsection add Add an algorithm
490
491 To add a new algorithm, one should check in the src/smpi/colls folder how other algorithms 
492 are coded. Using plain MPI code inside Simgrid can't be done, so algorithms have to be 
493 changed to use smpi version of the calls instead (MPI_Send will become smpi_mpi_send). Some functions may have different signatures than their MPI counterpart, please check the other algorithms or contact us using <a href="http://lists.gforge.inria.fr/mailman/listinfo/simgrid-devel">SimGrid developers mailing list</a>.
494
495 Example: adding a "pair" version of the Alltoall collective.
496
497  - Implement it in a file called alltoall-pair.c in the src/smpi/colls folder. This file should include colls_private.h.
498
499  - The name of the new algorithm function should be smpi_coll_tuned_alltoall_pair, with the same signature as MPI_Alltoall.
500
501  - Once the adaptation to SMPI code is done, add a reference to the file ("src/smpi/colls/alltoall-pair.c") in the SMPI_SRC part of the DefinePackages.cmake file inside buildtools/cmake, to allow the file to be built and distributed.
502
503  - To register the new version of the algorithm, simply add a line to the corresponding macro in src/smpi/colls/cools.h ( add a "COLL_APPLY(action, COLL_ALLTOALL_SIG, pair)" to the COLL_ALLTOALLS macro ). The algorithm should now be compiled and be selected when using --cfg=smpi/alltoall:pair at runtime.
504
505  - To add a test for the algorithm inside Simgrid's test suite, juste add the new algorithm name in the ALLTOALL_COLL list found inside buildtools/cmake/AddTests.cmake . When running ctest, a test for the new algorithm should be generated and executed. If it does not pass, please check your code or contact us.
506
507  - Feel free to push this new algorithm to the SMPI repository using Git.
508
509
510
511
512 */