Logo AND Algorithmique Numérique Distribuée

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