Logo AND Algorithmique Numérique Distribuée

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