Logo AND Algorithmique Numérique Distribuée

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