Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Merge remote-tracking branch 'upstream/master' into issue95
[simgrid.git] / docs / source / app_smpi.rst
1 .. _SMPI_doc:
2
3 ===============================
4 SMPI: Simulate MPI Applications
5 ===============================
6
7 .. raw:: html
8
9    <object id="TOC" data="graphical-toc.svg" type="image/svg+xml"></object>
10    <script>
11    window.onload=function() { // Wait for the SVG to be loaded before changing it
12      var elem=document.querySelector("#TOC").contentDocument.getElementById("SMPIBox")
13      elem.style="opacity:0.93999999;fill:#ff0000;fill-opacity:0.1";
14      elem.style="opacity:0.93999999;fill:#ff0000;fill-opacity:0.1;stroke:#000000;stroke-width:0.35277778;stroke-linecap:round;stroke-linejoin:round;stroke-miterlimit:4;stroke-dasharray:none;stroke-dashoffset:0;stroke-opacity:1";
15    }
16    </script>
17    <br/>
18    <br/>
19
20 SMPI enables the study of MPI application by emulating them on top of
21 the SimGrid simulator. This is particularly interesting to study
22 existing MPI applications within the comfort of the simulator.
23
24 To get started with SMPI, you should head to :ref:`the SMPI tutorial
25 <usecase_smpi>`. You may also want to read the `SMPI reference
26 article <https://hal.inria.fr/hal-01415484>`_ or these `introductory
27 slides <http://simgrid.org/tutorials/simgrid-smpi-101.pdf>`_.  If you
28 are new to MPI, you should first take our online `SMPI CourseWare
29 <https://simgrid.github.io/SMPI_CourseWare/>`_. It consists in several
30 projects that progressively introduce the MPI concepts. It proposes to
31 use SimGrid and SMPI to run the experiments, but the learning
32 objectives are centered on MPI itself.
33
34 Our goal is to enable the study of **unmodified MPI applications**.
35 Some constructs and features are still missing, but we can probably
36 add them on demand.  If you already used MPI before, SMPI should sound
37 very familiar to you: Use smpicc instead of mpicc, and smpirun instead
38 of mpirun. The main difference is that smpirun takes a :ref:`simulated
39 platform <platform>` as an extra parameter.
40
41 For **further scalability**, you may modify your code to speed up your
42 studies or save memory space.  Maximal **simulation accuracy**
43 requires some specific care from you.
44
45 .. _SMPI_online:
46
47 -----------------
48 Using SMPI online
49 -----------------
50
51 In this mode, your application is actually executed. Every computation
52 occurs for real while every communication is simulated. In addition,
53 the executions are automatically benchmarked so that their timings can
54 be applied within the simulator.
55
56 SMPI can also go offline by replaying a trace. :ref:`Trace replay
57 <SMPI_offline>` is usually ways faster than online simulation (because
58 the computation are skipped), but it can only applied to applications
59 with constant execution and communication patterns (for the exact same
60 reason).
61
62 ...................
63 Compiling your Code
64 ...................
65
66 If your application is in C, then simply use ``smpicc`` as a
67 compiler just like you use mpicc with other MPI implementations. This
68 script still calls your default compiler (gcc, clang, ...) and adds
69 the right compilation flags along the way. If your application is in
70 C++, Fortran 77 or Fortran 90, use respectively ``smpicxx``,
71 ``smpiff`` or ``smpif90``.
72
73 If you use cmake, set the variables ``MPI_C_COMPILER``, ``MPI_CXX_COMPILER`` and
74 ``MPI_Fortran_COMPILER`` to the full path of smpicc, smpicxx and smpiff (or
75 smpif90), respectively. Example:
76
77 .. code-block:: console
78
79    $ cmake -DMPI_C_COMPILER=/opt/simgrid/bin/smpicc -DMPI_CXX_COMPILER=/opt/simgrid/bin/smpicxx -DMPI_Fortran_COMPILER=/opt/simgrid/bin/smpiff .
80
81 ....................
82 Simulating your Code
83 ....................
84
85 Use the ``smpirun`` script as follows:
86
87 .. code-block:: console
88
89    $ smpirun -hostfile my_hostfile.txt -platform my_platform.xml ./program -blah
90
91 - ``my_hostfile.txt`` is a classical MPI hostfile (that is, this file
92   lists the machines on which the processes must be dispatched, one
93   per line)
94 - ``my_platform.xml`` is a classical SimGrid platform file. Of course,
95   the hosts of the hostfile must exist in the provided platform.
96 - ``./program`` is the MPI program to simulate, that you compiled with ``smpicc``
97 - ``-blah`` is a command-line parameter passed to this program.
98
99 ``smpirun`` accepts other parameters, such as ``-np`` if you don't
100 want to use all the hosts defined in the hostfile, ``-map`` to display
101 on which host each rank gets mapped of ``-trace`` to activate the
102 tracing during the simulation. You can get the full list by running
103 ``smpirun -help``
104
105 Finally, you can pass :ref:`any valid SimGrid parameter <options>` to your
106 program. In particular, you can pass ``--cfg=network/model:ns-3`` to
107 switch to use :ref:`model_ns3`. These parameters should be placed after
108 the name of your binary on the command line.
109
110 ...............................
111 Debugging your Code within SMPI
112 ...............................
113
114 If you want to explore the automatic platform and deployment files
115 that are generated by ``smpirun``, add ``-keep-temps`` to the command
116 line.
117
118 You can also run your simulation within valgrind or gdb using the
119 following commands. Once in GDB, each MPI ranks will be represented as
120 a regular thread, and you can explore the state of each of them as
121 usual.
122
123 .. code-block:: console
124
125    $ smpirun -wrapper valgrind ...other args...
126    $ smpirun -wrapper "gdb --args" --cfg=contexts/factory:thread ...other args...
127
128 Some shortcuts are available:
129
130 - ``-gdb`` is equivalent to ``-wrapper "gdb --args" -keep-temps``, to run within gdb debugger
131 - ``-lldb`` is equivalent to ``-wrapper "lldb --" -keep-temps``, to run within lldb debugger
132 - ``-vgdb`` is equivalent to ``-wrapper "valgrind --vgdb=yes --vgdb-error=0" -keep-temps``,
133   to run within valgrind and allow to attach a debugger
134
135 To help locate bottlenecks and largest allocations in the simulated application,
136 the -analyze flag can be passed to smpirun. It will activate
137 :ref:`smpi/display-timing<cfg=smpi/display-timing>` and
138 :ref:`smpi/display-allocs<cfg=smpi/display-allocs>` options and provide hints
139 at the end of execution.
140
141 SMPI will also report MPI handle (Comm, Request, Op, Datatype...) leaks
142 at the end of execution. This can help identify memory leaks that can trigger
143 crashes and slowdowns.
144 By default it only displays the number of leaked items detected.
145 Option :ref:`smpi/list-leaks:n<cfg=smpi/list-leaks>` can be used to display the
146 n first leaks encountered and their type. To get more information, running smpirun
147 with ``-wrapper "valgrind --leak-check=full --track-origins=yes"`` should show
148 the exact origin of leaked handles.
149 Known issue : MPI_Cancel may trigger internal leaks within SMPI.
150
151
152 .. _SMPI_use_colls:
153
154 ................................
155 Simulating Collective Operations
156 ................................
157
158 MPI collective operations are crucial to the performance of MPI
159 applications and must be carefully optimized according to many
160 parameters. Every existing implementation provides several algorithms
161 for each collective operation, and selects by default the best suited
162 one, depending on the sizes sent, the number of nodes, the
163 communicator, or the communication library being used.  These
164 decisions are based on empirical results and theoretical complexity
165 estimation, and are very different between MPI implementations. In
166 most cases, the users can also manually tune the algorithm used for
167 each collective operation.
168
169 SMPI can simulate the behavior of several MPI implementations:
170 OpenMPI, MPICH, `STAR-MPI <http://star-mpi.sourceforge.net/>`_, and
171 MVAPICH2. For that, it provides 115 collective algorithms and several
172 selector algorithms, that were collected directly in the source code
173 of the targeted MPI implementations.
174
175 You can switch the automatic selector through the
176 ``smpi/coll-selector`` configuration item. Possible values:
177
178  - **ompi:** default selection logic of OpenMPI (version 3.1.2)
179  - **mpich**: default selection logic of MPICH (version 3.3b)
180  - **mvapich2**: selection logic of MVAPICH2 (version 1.9) tuned
181    on the Stampede cluster
182  - **impi**: preliminary version of an Intel MPI selector (version
183    4.1.3, also tuned for the Stampede cluster). Due the closed source
184    nature of Intel MPI, some of the algorithms described in the
185    documentation are not available, and are replaced by mvapich ones.
186  - **default**: legacy algorithms used in the earlier days of
187    SimGrid. Do not use for serious perform performance studies.
188
189 .. todo:: default should not even exist.
190
191 ....................
192 Available Algorithms
193 ....................
194
195 You can also pick the algorithm used for each collective with the
196 corresponding configuration item. For example, to use the pairwise
197 alltoall algorithm, one should add ``--cfg=smpi/alltoall:pair`` to the
198 line. This will override the selector (if any) for this algorithm.  It
199 means that the selected algorithm will be used
200
201 .. Warning:: Some collective may require specific conditions to be
202    executed correctly (for instance having a communicator with a power
203    of two number of nodes only), which are currently not enforced by
204    Simgrid. Some crashes can be expected while trying these algorithms
205    with unusual sizes/parameters
206
207 MPI_Alltoall
208 ^^^^^^^^^^^^
209
210 Most of these are best described in `STAR-MPI's white paper <https://doi.org/10.1145/1183401.1183431>`_.
211
212 ``default``: naive one, by default. |br|
213 ``ompi``: use openmpi selector for the alltoall operations. |br|
214 ``mpich``: use mpich selector for the alltoall operations. |br|
215 ``mvapich2``: use mvapich2 selector for the alltoall operations. |br|
216 ``impi``: use intel mpi selector for the alltoall operations. |br|
217 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
218 ``bruck``: Described by Bruck et. al. in `this paper <http://ieeexplore.ieee.org/xpl/articleDetails.jsp?arnumber=642949>`_. |br|
219 ``2dmesh``: organizes the nodes as a two dimensional mesh, and perform allgather along the dimensions. |br|
220 ``3dmesh``: adds a third dimension to the previous algorithm. |br|
221 ``rdb``: recursive doubling``: extends the mesh to a nth dimension, each one containing two nodes. |br|
222 ``pair``: pairwise exchange, only works for power of 2 procs, size-1 steps, each process sends and receives from the same process at each step. |br|
223 ``pair_light_barrier``: same, with small barriers between steps to avoid contention. |br|
224 ``pair_mpi_barrier``: same, with MPI_Barrier used. |br|
225 ``pair_one_barrier``: only one barrier at the beginning. |br|
226 ``ring``: size-1 steps, at each step a process send to process (n+i)%size, and receives from (n-i)%size. |br|
227 ``ring_light_barrier``: same, with small barriers between some phases to avoid contention. |br|
228 ``ring_mpi_barrier``: same, with MPI_Barrier used. |br|
229 ``ring_one_barrier``: only one barrier at the beginning. |br|
230 ``basic_linear``: posts all receives and all sends, starts the communications, and waits for all communication to finish. |br|
231 ``mvapich2_scatter_dest``: isend/irecv with scattered destinations, posting only a few messages at the same time. |br|
232
233 MPI_Alltoallv
234 ^^^^^^^^^^^^^
235
236 ``default``: naive one, by default. |br|
237 ``ompi``: use openmpi selector for the alltoallv operations. |br|
238 ``mpich``: use mpich selector for the alltoallv operations. |br|
239 ``mvapich2``: use mvapich2 selector for the alltoallv operations. |br|
240 ``impi``: use intel mpi selector for the alltoallv operations. |br|
241 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
242 ``bruck``: same as alltoall. |br|
243 ``pair``: same as alltoall. |br|
244 ``pair_light_barrier``: same as alltoall. |br|
245 ``pair_mpi_barrier``: same as alltoall. |br|
246 ``pair_one_barrier``: same as alltoall. |br|
247 ``ring``: same as alltoall. |br|
248 ``ring_light_barrier``: same as alltoall. |br|
249 ``ring_mpi_barrier``: same as alltoall. |br|
250 ``ring_one_barrier``: same as alltoall. |br|
251 ``ompi_basic_linear``: same as alltoall. |br|
252
253 MPI_Gather
254 ^^^^^^^^^^
255
256 ``default``: naive one, by default. |br|
257 ``ompi``: use openmpi selector for the gather operations. |br|
258 ``mpich``: use mpich selector for the gather operations. |br|
259 ``mvapich2``: use mvapich2 selector for the gather operations. |br|
260 ``impi``: use intel mpi selector for the gather operations. |br|
261 ``automatic (experimental)``: use an automatic self-benchmarking algorithm which will iterate over all implemented versions and output the best. |br|
262 ``ompi_basic_linear``: basic linear algorithm from openmpi, each process sends to the root. |br|
263 ``ompi_binomial``: binomial tree algorithm. |br|
264 ``ompi_linear_sync``: same as basic linear, but with a synchronization at the beginning and message cut into two segments. |br|
265 ``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. |br|
266
267 MPI_Barrier
268 ^^^^^^^^^^^
269
270 ``default``: naive one, by default. |br|
271 ``ompi``: use openmpi selector for the barrier operations. |br|
272 ``mpich``: use mpich selector for the barrier operations. |br|
273 ``mvapich2``: use mvapich2 selector for the barrier operations. |br|
274 ``impi``: use intel mpi selector for the barrier operations. |br|
275 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
276 ``ompi_basic_linear``: all processes send to root. |br|
277 ``ompi_two_procs``: special case for two processes. |br|
278 ``ompi_bruck``: nsteps = sqrt(size), at each step, exchange data with rank-2^k and rank+2^k. |br|
279 ``ompi_recursivedoubling``: recursive doubling algorithm. |br|
280 ``ompi_tree``: recursive doubling type algorithm, with tree structure. |br|
281 ``ompi_doublering``: double ring algorithm. |br|
282 ``mvapich2_pair``: pairwise algorithm. |br|
283 ``mpich_smp``: barrier intra-node, then inter-node. |br|
284
285 MPI_Scatter
286 ^^^^^^^^^^^
287
288 ``default``: naive one, by default. |br|
289 ``ompi``: use openmpi selector for the scatter operations. |br|
290 ``mpich``: use mpich selector for the scatter operations. |br|
291 ``mvapich2``: use mvapich2 selector for the scatter operations. |br|
292 ``impi``: use intel mpi selector for the scatter operations. |br|
293 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
294 ``ompi_basic_linear``: basic linear scatter. |br|
295 ``ompi_binomial``: binomial tree scatter. |br|
296 ``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. |br|
297 ``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. |br|
298
299 MPI_Reduce
300 ^^^^^^^^^^
301
302 ``default``: naive one, by default. |br|
303 ``ompi``: use openmpi selector for the reduce operations. |br|
304 ``mpich``: use mpich selector for the reduce operations. |br|
305 ``mvapich2``: use mvapich2 selector for the reduce operations. |br|
306 ``impi``: use intel mpi selector for the reduce operations. |br|
307 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
308 ``arrival_pattern_aware``: root exchanges with the first process to arrive. |br|
309 ``binomial``: uses a binomial tree. |br|
310 ``flat_tree``: uses a flat tree. |br|
311 ``NTSL``: Non-topology-specific pipelined linear-bcast function. |br| 0->1, 1->2 ,2->3, ....., ->last node: in a pipeline fashion, with segments of 8192 bytes. |br|
312 ``scatter_gather``: scatter then gather. |br|
313 ``ompi_chain``: openmpi reduce algorithms are built on the same basis, but the topology is generated differently for each flavor. chain = chain with spacing of size/2, and segment size of 64KB. |br|
314 ``ompi_pipeline``: same with pipeline (chain with spacing of 1), segment size depends on the communicator size and the message size. |br|
315 ``ompi_binary``: same with binary tree, segment size of 32KB. |br|
316 ``ompi_in_order_binary``: same with binary tree, enforcing order on the operations. |br|
317 ``ompi_binomial``: same with binomial algo (redundant with default binomial one in most cases). |br|
318 ``ompi_basic_linear``: basic algorithm, each process sends to root. |br|
319 ``mvapich2_knomial``: k-nomial algorithm. Default factor is 4 (mvapich2 selector adapts it through tuning). |br|
320 ``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. |br|
321 ``rab``: `Rabenseifner <https://fs.hlrs.de/projects/par/mpi//myreduce.html>`_'s reduce algorithm. |br|
322
323 MPI_Allreduce
324 ^^^^^^^^^^^^^
325
326 ``default``: naive one, by defautl. |br|
327 ``ompi``: use openmpi selector for the allreduce operations. |br|
328 ``mpich``: use mpich selector for the allreduce operations. |br|
329 ``mvapich2``: use mvapich2 selector for the allreduce operations. |br|
330 ``impi``: use intel mpi selector for the allreduce operations. |br|
331 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
332 ``lr``: logical ring reduce-scatter then logical ring allgather. |br|
333 ``rab1``: variations of the  `Rabenseifner <https://fs.hlrs.de/projects/par/mpi//myreduce.html>`_ algorithm: reduce_scatter then allgather. |br|
334 ``rab2``: variations of the  `Rabenseifner <https://fs.hlrs.de/projects/par/mpi//myreduce.html>`_ algorithm: alltoall then allgather. |br|
335 ``rab_rsag``: variation of the  `Rabenseifner <https://fs.hlrs.de/projects/par/mpi//myreduce.html>`_ algorithm: recursive doubling reduce_scatter then recursive doubling allgather. |br|
336 ``rdb``: recursive doubling. |br|
337 ``smp_binomial``: binomial tree with smp: binomial intra. |br| SMP reduce, inter reduce, inter broadcast then intra broadcast. |br|
338 ``smp_binomial_pipeline``: same with segment size = 4096 bytes. |br|
339 ``smp_rdb``: intra``: binomial allreduce, inter: Recursive doubling allreduce, intra``: binomial broadcast. |br|
340 ``smp_rsag``: intra: binomial allreduce, inter: reduce-scatter, inter:allgather, intra: binomial broadcast. |br|
341 ``smp_rsag_lr``: intra: binomial allreduce, inter: logical ring reduce-scatter, logical ring inter:allgather, intra: binomial broadcast. |br|
342 ``smp_rsag_rab``: intra: binomial allreduce, inter: rab reduce-scatter, rab inter:allgather, intra: binomial broadcast. |br|
343 ``redbcast``: reduce then broadcast, using default or tuned algorithms if specified. |br|
344 ``ompi_ring_segmented``: ring algorithm used by OpenMPI. |br|
345 ``mvapich2_rs``: rdb for small messages, reduce-scatter then allgather else. |br|
346 ``mvapich2_two_level``: SMP-aware algorithm, with mpich as intra algorithm, and rdb as inter (Change this behavior by using mvapich2 selector to use tuned values). |br|
347 ``rab``: default `Rabenseifner <https://fs.hlrs.de/projects/par/mpi//myreduce.html>`_ implementation. |br|
348
349 MPI_Reduce_scatter
350 ^^^^^^^^^^^^^^^^^^
351
352 ``default``: naive one, by default. |br|
353 ``ompi``: use openmpi selector for the reduce_scatter operations. |br|
354 ``mpich``: use mpich selector for the reduce_scatter operations. |br|
355 ``mvapich2``: use mvapich2 selector for the reduce_scatter operations. |br|
356 ``impi``: use intel mpi selector for the reduce_scatter operations. |br|
357 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
358 ``ompi_basic_recursivehalving``: recursive halving version from OpenMPI. |br|
359 ``ompi_ring``: ring version from OpenMPI. |br|
360 ``mpich_pair``: pairwise exchange version from MPICH. |br|
361 ``mpich_rdb``: recursive doubling version from MPICH. |br|
362 ``mpich_noncomm``: only works for power of 2 procs, recursive doubling for noncommutative ops. |br|
363
364
365 MPI_Allgather
366 ^^^^^^^^^^^^^
367
368 ``default``: naive one, by default. |br|
369 ``ompi``: use openmpi selector for the allgather operations. |br|
370 ``mpich``: use mpich selector for the allgather operations. |br|
371 ``mvapich2``: use mvapich2 selector for the allgather operations. |br|
372 ``impi``: use intel mpi selector for the allgather operations. |br|
373 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
374 ``2dmesh``: see alltoall. |br|
375 ``3dmesh``: see alltoall. |br|
376 ``bruck``: Described by Bruck et.al. in <a href="http://ieeexplore.ieee.org/xpl/articleDetails.jsp?arnumber=642949"> Efficient algorithms for all-to-all communications in multiport message-passing systems</a>. |br|
377 ``GB``: Gather - Broadcast (uses tuned version if specified). |br|
378 ``loosely_lr``: Logical Ring with grouping by core (hardcoded, default processes/node: 4). |br|
379 ``NTSLR``: Non Topology Specific Logical Ring. |br|
380 ``NTSLR_NB``: Non Topology Specific Logical Ring, Non Blocking operations. |br|
381 ``pair``: see alltoall. |br|
382 ``rdb``: see alltoall. |br|
383 ``rhv``: only power of 2 number of processes. |br|
384 ``ring``: see alltoall. |br|
385 ``SMP_NTS``: gather to root of each SMP, then every root of each SMP node. post INTER-SMP Sendrecv, then do INTRA-SMP Bcast for each receiving message, using logical ring algorithm (hardcoded, default processes/SMP: 8). |br|
386 ``smp_simple``: gather to root of each SMP, then every root of each SMP node post INTER-SMP Sendrecv, then do INTRA-SMP Bcast for each receiving message, using simple algorithm (hardcoded, default processes/SMP: 8). |br|
387 ``spreading_simple``: from node i, order of communications is i -> i + 1, i -> i + 2, ..., i -> (i + p -1) % P. |br|
388 ``ompi_neighborexchange``: Neighbor Exchange algorithm for allgather. Described by Chen et.al. in  `Performance Evaluation of Allgather Algorithms on Terascale Linux Cluster with Fast Ethernet <http://ieeexplore.ieee.org/xpl/articleDetails.jsp?tp=&arnumber=1592302>`_. |br|
389 ``mvapich2_smp``: SMP aware algorithm, performing intra-node gather, inter-node allgather with one process/node, and bcast intra-node
390
391 MPI_Allgatherv
392 ^^^^^^^^^^^^^^
393
394 ``default``: naive one, by default. |br|
395 ``ompi``: use openmpi selector for the allgatherv operations. |br|
396 ``mpich``: use mpich selector for the allgatherv operations. |br|
397 ``mvapich2``: use mvapich2 selector for the allgatherv operations. |br|
398 ``impi``: use intel mpi selector for the allgatherv operations. |br|
399 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
400 ``GB``: Gatherv - Broadcast (uses tuned version if specified, but only for Bcast, gatherv is not tuned). |br|
401 ``pair``: see alltoall. |br|
402 ``ring``: see alltoall. |br|
403 ``ompi_neighborexchange``: see allgather. |br|
404 ``ompi_bruck``: see allgather. |br|
405 ``mpich_rdb``: recursive doubling algorithm from MPICH. |br|
406 ``mpich_ring``: ring algorithm from MPICh - performs differently from the  one from STAR-MPI.
407
408 MPI_Bcast
409 ^^^^^^^^^
410
411 ``default``: naive one, by default. |br|
412 ``ompi``: use openmpi selector for the bcast operations. |br|
413 ``mpich``: use mpich selector for the bcast operations. |br|
414 ``mvapich2``: use mvapich2 selector for the bcast operations. |br|
415 ``impi``: use intel mpi selector for the bcast operations. |br|
416 ``automatic (experimental)``: use an automatic self-benchmarking algorithm. |br|
417 ``arrival_pattern_aware``: root exchanges with the first process to arrive. |br|
418 ``arrival_pattern_aware_wait``: same with slight variation. |br|
419 ``binomial_tree``: binomial tree exchange. |br|
420 ``flattree``: flat tree exchange. |br|
421 ``flattree_pipeline``: flat tree exchange, message split into 8192 bytes pieces. |br|
422 ``NTSB``: Non-topology-specific pipelined binary tree with 8192 bytes pieces. |br|
423 ``NTSL``: Non-topology-specific pipelined linear with 8192 bytes pieces. |br|
424 ``NTSL_Isend``: Non-topology-specific pipelined linear with 8192 bytes pieces, asynchronous communications. |br|
425 ``scatter_LR_allgather``: scatter followed by logical ring allgather. |br|
426 ``scatter_rdb_allgather``: scatter followed by recursive doubling allgather. |br|
427 ``arrival_scatter``: arrival pattern aware scatter-allgather. |br|
428 ``SMP_binary``: binary tree algorithm with 8 cores/SMP. |br|
429 ``SMP_binomial``: binomial tree algorithm with 8 cores/SMP. |br|
430 ``SMP_linear``: linear algorithm with 8 cores/SMP. |br|
431 ``ompi_split_bintree``: binary tree algorithm from OpenMPI, with message split in 8192 bytes pieces. |br|
432 ``ompi_pipeline``: pipeline algorithm from OpenMPI, with message split in 128KB pieces. |br|
433 ``mvapich2_inter_node``: Inter node default mvapich worker. |br|
434 ``mvapich2_intra_node``: Intra node default mvapich worker. |br|
435 ``mvapich2_knomial_intra_node``:  k-nomial intra node default mvapich worker. default factor is 4.
436
437 Automatic Evaluation
438 ^^^^^^^^^^^^^^^^^^^^
439
440 .. warning:: This is still very experimental.
441
442 An automatic version is available for each collective (or even as a selector). This specific
443 version will loop over all other implemented algorithm for this particular collective, and apply
444 them while benchmarking the time taken for each process. It will then output the quickest for
445 each process, and the global quickest. This is still unstable, and a few algorithms which need
446 specific number of nodes may crash.
447
448 Adding an algorithm
449 ^^^^^^^^^^^^^^^^^^^
450
451 To add a new algorithm, one should check in the src/smpi/colls folder
452 how other algorithms are coded. Using plain MPI code inside Simgrid
453 can't be done, so algorithms have to be changed to use smpi version of
454 the calls instead (MPI_Send will become smpi_mpi_send). Some functions
455 may have different signatures than their MPI counterpart, please check
456 the other algorithms or contact us using the `>SimGrid
457 developers mailing list <http://lists.gforge.inria.fr/mailman/listinfo/simgrid-devel>`_.
458
459 Example: adding a "pair" version of the Alltoall collective.
460
461  - Implement it in a file called alltoall-pair.c in the src/smpi/colls folder. This file should include colls_private.hpp.
462
463  - The name of the new algorithm function should be smpi_coll_tuned_alltoall_pair, with the same signature as MPI_Alltoall.
464
465  - 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.
466
467  - 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.
468
469  - 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 teshsuite/smpi/CMakeLists.txt . 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.
470
471  - Please submit your patch for inclusion in SMPI, for example through a pull request on GitHub or directly per email.
472
473
474 Tracing of Internal Communications
475 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
476
477 By default, the collective operations are traced as a unique operation
478 because tracing all point-to-point communications composing them could
479 result in overloaded, hard to interpret traces. If you want to debug
480 and compare collective algorithms, you should set the
481 ``tracing/smpi/internals`` configuration item to 1 instead of 0.
482
483 Here are examples of two alltoall collective algorithms runs on 16 nodes,
484 the first one with a ring algorithm, the second with a pairwise one.
485
486 .. image:: /img/smpi_simgrid_alltoall_ring_16.png
487    :align: center
488
489 Alltoall on 16 Nodes with the Ring Algorithm.
490
491 .. image:: /img/smpi_simgrid_alltoall_pair_16.png
492    :align: center
493
494 Alltoall on 16 Nodes with the Pairwise Algorithm.
495
496 -------------------------
497 What can run within SMPI?
498 -------------------------
499
500 You can run unmodified MPI applications (both C/C++ and Fortran) within
501 SMPI, provided that you only use MPI calls that we implemented. Global
502 variables should be handled correctly on Linux systems.
503
504 ....................
505 MPI coverage of SMPI
506 ....................
507
508 SMPI support a large faction of the MPI interface: we pass many of the MPICH coverage tests, and many of the existing
509 :ref:`proxy apps <SMPI_proxy_apps>` run almost unmodified on top of SMPI. But our support is still incomplete, with I/O
510 primitives the being one of the major missing feature.
511
512 The full list of functions that remain to be implemented is documented in the file `include/smpi/smpi.h
513 <https://framagit.org/simgrid/simgrid/tree/master/include/smpi/smpi.h>`_ in your version of SimGrid, between two lines
514 containing the ``FIXME`` marker. If you miss a feature, please get in touch with us: we can guide you through the SimGrid
515 code to help you implementing it, and we'd be glad to integrate your contribution to the main project.
516
517 .. _SMPI_what_globals:
518
519 .................................
520 Privatization of global variables
521 .................................
522
523 Concerning the globals, the problem comes from the fact that usually,
524 MPI processes run as real UNIX processes while they are all folded
525 into threads of a unique system process in SMPI. Global variables are
526 usually private to each MPI process while they become shared between
527 the processes in SMPI.  The problem and some potential solutions are
528 discussed in this article: `Automatic Handling of Global Variables for
529 Multi-threaded MPI Programs
530 <http://charm.cs.illinois.edu/newPapers/11-23/paper.pdf>` (note that
531 this article does not deal with SMPI but with a competing solution
532 called AMPI that suffers of the same issue).  This point used to be
533 problematic in SimGrid, but the problem should now be handled
534 automatically on Linux.
535
536 Older versions of SimGrid came with a script that automatically
537 privatized the globals through static analysis of the source code. But
538 our implementation was not robust enough to be used in production, so
539 it was removed at some point. Currently, SMPI comes with two
540 privatization mechanisms that you can :ref:`select at runtime
541 <cfg=smpi/privatization>`.  The dlopen approach is used by
542 default as it is much faster and still very robust.  The mmap approach
543 is an older approach that proves to be slower.
544
545 With the **mmap approach**, SMPI duplicates and dynamically switch the
546 ``.data`` and ``.bss`` segments of the ELF process when switching the
547 MPI ranks. This allows each ranks to have its own copy of the global
548 variables.  No copy actually occurs as this mechanism uses ``mmap()``
549 for efficiency. This mechanism is considered to be very robust on all
550 systems supporting ``mmap()`` (Linux and most BSDs). Its performance
551 is questionable since each context switch between MPI ranks induces
552 several syscalls to change the ``mmap`` that redirects the ``.data``
553 and ``.bss`` segments to the copies of the new rank. The code will
554 also be copied several times in memory, inducing a slight increase of
555 memory occupation.
556
557 Another limitation is that SMPI only accounts for global variables
558 defined in the executable. If the processes use external global
559 variables from dynamic libraries, they won't be switched
560 correctly. The easiest way to solve this is to statically link against
561 the library with these globals. This way, each MPI rank will get its
562 own copy of these libraries. Of course you should never statically
563 link against the SimGrid library itself.
564
565 With the **dlopen approach**, SMPI loads several copies of the same
566 executable in memory as if it were a library, so that the global
567 variables get naturally duplicated. It first requires the executable
568 to be compiled as a relocatable binary, which is less common for
569 programs than for libraries. But most distributions are now compiled
570 this way for security reason as it allows one to randomize the address
571 space layout. It should thus be safe to compile most (any?) program
572 this way.  The second trick is that the dynamic linker refuses to link
573 the exact same file several times, be it a library or a relocatable
574 executable. It makes perfectly sense in the general case, but we need
575 to circumvent this rule of thumb in our case. To that extend, the
576 binary is copied in a temporary file before being re-linked against.
577 ``dlmopen()`` cannot be used as it only allows 256 contextes, and as it
578 would also duplicate simgrid itself.
579
580 This approach greatly speeds up the context switching, down to about
581 40 CPU cycles with our raw contextes, instead of requesting several
582 syscalls with the ``mmap()`` approach. Another advantage is that it
583 permits one to run the SMPI contexts in parallel, which is obviously not
584 possible with the ``mmap()`` approach. It was tricky to implement, but
585 we are not aware of any flaws, so smpirun activates it by default.
586
587 In the future, it may be possible to further reduce the memory and
588 disk consumption. It seems that we could `punch holes
589 <https://lwn.net/Articles/415889/>`_ in the files before dl-loading
590 them to remove the code and constants, and mmap these area onto a
591 unique copy. If done correctly, this would reduce the disk- and
592 memory- usage to the bare minimum, and would also reduce the pressure
593 on the CPU instruction cache. See the `relevant bug
594 <https://github.com/simgrid/simgrid/issues/137>`_ on github for
595 implementation leads.\n
596
597 Also, currently, only the binary is copied and dlopen-ed for each MPI
598 rank. We could probably extend this to external dependencies, but for
599 now, any external dependencies must be statically linked into your
600 application. As usual, simgrid itself shall never be statically linked
601 in your app. You don't want to give a copy of SimGrid to each MPI rank:
602 that's ways too much for them to deal with.
603
604 .. todo: speak of smpi/privatize-libs here
605
606 ----------------------------------------------
607 Adapting your MPI code for further scalability
608 ----------------------------------------------
609
610 As detailed in the `reference article
611 <http://hal.inria.fr/hal-01415484>`_, you may want to adapt your code
612 to improve the simulation performance. But these tricks may seriously
613 hinder the result quality (or even prevent the app to run) if used
614 wrongly. We assume that if you want to simulate an HPC application,
615 you know what you are doing. Don't prove us wrong!
616
617 ..............................
618 Reducing your memory footprint
619 ..............................
620
621 If you get short on memory (the whole app is executed on a single node when
622 simulated), you should have a look at the SMPI_SHARED_MALLOC and
623 SMPI_SHARED_FREE macros. It allows one to share memory areas between processes: The
624 purpose of these macro is that the same line malloc on each process will point
625 to the exact same memory area. So if you have a malloc of 2M and you have 16
626 processes, this macro will change your memory consumption from 2M*16 to 2M
627 only. Only one block for all processes.
628
629 If your program is ok with a block containing garbage value because all
630 processes write and read to the same place without any kind of coordination,
631 then this macro can dramatically shrink your memory consumption. For example,
632 that will be very beneficial to a matrix multiplication code, as all blocks will
633 be stored on the same area. Of course, the resulting computations will useless,
634 but you can still study the application behavior this way.
635
636 Naturally, this won't work if your code is data-dependent. For example, a Jacobi
637 iterative computation depends on the result computed by the code to detect
638 convergence conditions, so turning them into garbage by sharing the same memory
639 area between processes does not seem very wise. You cannot use the
640 SMPI_SHARED_MALLOC macro in this case, sorry.
641
642 This feature is demoed by the example file
643 `examples/smpi/NAS/dt.c <https://framagit.org/simgrid/simgrid/tree/master/examples/smpi/NAS/dt.c>`_
644
645 .. _SMPI_use_faster:
646
647 .........................
648 Toward Faster Simulations
649 .........................
650
651 If your application is too slow, try using SMPI_SAMPLE_LOCAL,
652 SMPI_SAMPLE_GLOBAL and friends to indicate which computation loops can
653 be sampled. Some of the loop iterations will be executed to measure
654 their duration, and this duration will be used for the subsequent
655 iterations. These samples are done per processor with
656 SMPI_SAMPLE_LOCAL, and shared between all processors with
657 SMPI_SAMPLE_GLOBAL. Of course, none of this will work if the execution
658 time of your loop iteration are not stable. If some parameters have an
659 incidence on the timing of a kernel, and if they are reused often
660 (same kernel launched with a few different sizes during the run, for example),
661 SMPI_SAMPLE_LOCAL_TAG and SMPI_SAMPLE_GLOBAL_TAG can be used, with a tag
662 as last parameter, to differentiate between calls. The tag is a character
663 chain crafted by the user, with a maximum size of 128, and should include
664 what is necessary to group calls of a given size together.
665
666 This feature is demoed by the example file
667 `examples/smpi/NAS/ep.c <https://framagit.org/simgrid/simgrid/tree/master/examples/smpi/NAS/ep.c>`_
668
669 .............................
670 Ensuring Accurate Simulations
671 .............................
672
673 Out of the box, SimGrid may give you fairly accurate results, but
674 there is a plenty of factors that could go wrong and make your results
675 inaccurate or even plainly wrong. Actually, you can only get accurate
676 results of a nicely built model, including both the system hardware
677 and your application. Such models are hard to pass over and reuse in
678 other settings, because elements that are not relevant to an
679 application (say, the latency of point-to-point communications,
680 collective operation implementation details or CPU-network
681 interaction) may be irrelevant to another application. The dream of
682 the perfect model, encompassing every aspects is only a chimera, as
683 the only perfect model of the reality is the reality. If you go for
684 simulation, then you have to ignore some irrelevant aspects of the
685 reality, but which aspects are irrelevant is actually
686 application-dependent...
687
688 The only way to assess whether your settings provide accurate results
689 is to double-check these results. If possible, you should first run
690 the same experiment in simulation and in real life, gathering as much
691 information as you can. Try to understand the discrepancies in the
692 results that you observe between both settings (visualization can be
693 precious for that). Then, try to modify your model (of the platform,
694 of the collective operations) to reduce the most preeminent differences.
695
696 If the discrepancies come from the computing time, try adapting the
697 ``smpi/host-speed``: reduce it if your simulation runs faster than in
698 reality. If the error come from the communication, then you need to
699 fiddle with your platform file.
700
701 Be inventive in your modeling. Don't be afraid if the names given by
702 SimGrid does not match the real names: we got very good results by
703 modeling multicore/GPU machines with a set of separate hosts
704 interconnected with very fast networks (but don't trust your model
705 because it has the right names in the right place either).
706
707 Finally, you may want to check `this article
708 <https://hal.inria.fr/hal-00907887>`_ on the classical pitfalls in
709 modeling distributed systems.
710
711 .. _SMPI_proxy_apps:
712
713 ----------------------
714 Examples of SMPI Usage
715 ----------------------
716
717 A small amount of examples can be found directly in the SimGrid
718 archive, under `examples/smpi <https://framagit.org/simgrid/simgrid/-/tree/master/examples/smpi>`_.
719 Some show how to simply run MPI code in SimGrid, how to use the
720 tracing/replay mechanism or how to use plugins written in S4U to
721 extend the simulator abilities.
722
723 Another source of examples lay in the SimGrid archive, under
724 `teshsuite/smpi <https://framagit.org/simgrid/simgrid/-/tree/master/examples/smpi>`_.
725 They are not in the ``examples`` directory because they probably don't
726 constitute pedagogical examples. Instead, they are intended to stress
727 our implementation during the tests. Some of you may be interested
728 anyway.
729
730 But the best source of SMPI examples is certainly the `proxy app
731 <https://framagit.org/simgrid/SMPI-proxy-apps>`_ external project.
732 Proxy apps are scale models of real, massive HPC applications: each of
733 them exhibits the same communication and computation patterns than the
734 massive application that it stands for. But they last only a few
735 thousands lines instead of some millions of lines. These proxy apps
736 are usually provided for educational purpose, and also to ensure that
737 the represented large HPC applications will correctly work with the
738 next generation of runtimes and hardware. `This project
739 <https://framagit.org/simgrid/SMPI-proxy-apps>`_ gathers proxy apps
740 from different sources, along with the patches needed (if any) to run
741 them on top of SMPI.
742
743 -------------------------
744 Troubleshooting with SMPI
745 -------------------------
746
747 .................................
748 ./configure refuses to use smpicc
749 .................................
750
751 If your ``./configure`` reports that the compiler is not
752 functional or that you are cross-compiling, try to define the
753 ``SMPI_PRETEND_CC`` environment variable before running the
754 configuration.
755
756 .. code-block:: console
757
758    $ SMPI_PRETEND_CC=1 ./configure # here come the configure parameters
759    $ make
760
761 Indeed, the programs compiled with ``smpicc`` cannot be executed
762 without ``smpirun`` (they are shared libraries and do weird things on
763 startup), while configure wants to test them directly.  With
764 ``SMPI_PRETEND_CC`` smpicc does not compile as shared, and the SMPI
765 initialization stops and returns 0 before doing anything that would
766 fail without ``smpirun``.
767
768 .. warning::
769
770   Make sure that SMPI_PRETEND_CC is only set when calling ./configure,
771   not during the actual execution, or any program compiled with smpicc
772   will stop before starting.
773
774 ..............................................
775 ./configure does not pick smpicc as a compiler
776 ..............................................
777
778 In addition to the previous answers, some projects also need to be
779 explicitly told what compiler to use, as follows:
780
781 .. code-block:: console
782
783    $ SMPI_PRETEND_CC=1 ./configure CC=smpicc # here come the other configure parameters
784    $ make
785
786 Maybe your configure is using another variable, such as ``cc`` (in
787 lower case) or similar. Just check the logs.
788
789 .....................................
790 error: unknown type name 'useconds_t'
791 .....................................
792
793 Try to add ``-D_GNU_SOURCE`` to your compilation line to get rid
794 of that error.
795
796 The reason is that SMPI provides its own version of ``usleep(3)``
797 to override it and to block in the simulation world, not in the real
798 one. It needs the ``useconds_t`` type for that, which is declared
799 only if you declare ``_GNU_SOURCE`` before including
800 ``unistd.h``. If your project includes that header file before
801 SMPI, then you need to ensure that you pass the right configuration
802 defines as advised above.
803
804
805
806 .. _SMPI_offline:
807
808 -----------------------------
809 Trace Replay and Offline SMPI
810 -----------------------------
811
812 Although SMPI is often used for :ref:`online simulation
813 <SMPI_online>`, where the application is executed for real, you can
814 also go for offline simulation through trace replay.
815
816 SimGrid uses time-independent traces, in which each actor is given a
817 script of the actions to do sequentially. These trace files can
818 actually be captured with the online version of SMPI, as follows:
819
820 .. code-block:: console
821
822    $ smpirun -trace-ti --cfg=tracing/filename:LU.A.32 -np 32 -platform ../cluster_backbone.xml bin/lu.A.32
823
824 The produced trace is composed of a file ``LU.A.32`` and a folder
825 ``LU.A.32_files``. The file names don't match with the MPI ranks, but
826 that's expected.
827
828 To replay this with SMPI, you need to first compile the provided
829 ``smpi_replay.cpp`` file, that comes from
830 `simgrid/examples/smpi/replay
831 <https://framagit.org/simgrid/simgrid/tree/master/examples/smpi/replay>`_.
832
833 .. code-block:: console
834
835    $ smpicxx ../replay.cpp -O3 -o ../smpi_replay
836
837 Afterward, you can replay your trace in SMPI as follows:
838
839 .. code-block:: console
840
841    $ smpirun -np 32 -platform ../cluster_torus.xml -ext smpi_replay ../smpi_replay LU.A.32
842
843 All the outputs are gone, as the application is not really simulated
844 here. Its trace is simply replayed. But if you visualize the live
845 simulation and the replay, you will see that the behavior is
846 unchanged. The simulation does not run much faster on this very
847 example, but this becomes very interesting when your application
848 is computationally hungry.
849
850 .. |br| raw:: html
851
852    <br />